xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision 940cd9d6a860b4be8b5cd21e391a8cea58b4a45b)
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 
45*940cd9d6SSatish Balay /* declare MumpsScalar */
46*940cd9d6SSatish Balay #if defined(PETSC_USE_COMPLEX)
47*940cd9d6SSatish Balay #if defined(PETSC_USE_REAL_SINGLE)
48*940cd9d6SSatish Balay #define MumpsScalar mumps_complex
49*940cd9d6SSatish Balay #else
50*940cd9d6SSatish Balay #define MumpsScalar mumps_double_complex
51*940cd9d6SSatish Balay #endif
52*940cd9d6SSatish Balay #else
53*940cd9d6SSatish Balay #define MumpsScalar PetscScalar
54*940cd9d6SSatish Balay #endif
553d472b54SHong Zhang 
56397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */
57397b6df1SKris Buschelman #define ICNTL(I) icntl[(I)-1]
58397b6df1SKris Buschelman #define CNTL(I) cntl[(I)-1]
59397b6df1SKris Buschelman #define INFOG(I) infog[(I)-1]
60a7aca84bSHong Zhang #define INFO(I) info[(I)-1]
61397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1]
62adc1d99fSHong Zhang #define RINFO(I) rinfo[(I)-1]
63397b6df1SKris Buschelman 
64397b6df1SKris Buschelman typedef struct {
65397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
662907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
672907cef9SHong Zhang   CMUMPS_STRUC_C id;
682907cef9SHong Zhang #else
69397b6df1SKris Buschelman   ZMUMPS_STRUC_C id;
702907cef9SHong Zhang #endif
712907cef9SHong Zhang #else
722907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
732907cef9SHong Zhang   SMUMPS_STRUC_C id;
74397b6df1SKris Buschelman #else
75397b6df1SKris Buschelman   DMUMPS_STRUC_C id;
76397b6df1SKris Buschelman #endif
772907cef9SHong Zhang #endif
782907cef9SHong Zhang 
79397b6df1SKris Buschelman   MatStructure matstruc;
80c1490034SHong Zhang   PetscMPIInt  myid,size;
81a5e57a09SHong Zhang   PetscInt     *irn,*jcn,nz,sym;
82397b6df1SKris Buschelman   PetscScalar  *val;
83397b6df1SKris Buschelman   MPI_Comm     comm_mumps;
8464e6c443SBarry Smith   PetscBool    isAIJ,CleanUpMUMPS;
85a5e57a09SHong Zhang   PetscInt     ICNTL9_pre;           /* check if ICNTL(9) is changed from previous MatSolve */
86801fbe65SHong Zhang   VecScatter   scat_rhs, scat_sol;   /* used by MatSolve() */
87801fbe65SHong Zhang   Vec          b_seq,x_seq;
88b34f08ffSHong Zhang   PetscInt     ninfo,*info;          /* display INFO */
892205254eSKarl Rupp 
90bf0cc555SLisandro Dalcin   PetscErrorCode (*Destroy)(Mat);
91bccb9932SShri Abhyankar   PetscErrorCode (*ConvertToTriples)(Mat, int, MatReuse, int*, int**, int**, PetscScalar**);
92f0c56d0fSKris Buschelman } Mat_MUMPS;
93f0c56d0fSKris Buschelman 
9409573ac7SBarry Smith extern PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);
95b24902e0SBarry Smith 
96397b6df1SKris Buschelman /*
97d341cd04SHong Zhang   MatConvertToTriples_A_B - convert Petsc matrix to triples: row[nz], col[nz], val[nz]
98d341cd04SHong Zhang 
99397b6df1SKris Buschelman   input:
10067877ebaSShri Abhyankar     A       - matrix in aij,baij or sbaij (bs=1) format
101397b6df1SKris Buschelman     shift   - 0: C style output triple; 1: Fortran style output triple.
102bccb9932SShri Abhyankar     reuse   - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple
103bccb9932SShri Abhyankar               MAT_REUSE_MATRIX:   only the values in v array are updated
104397b6df1SKris Buschelman   output:
105397b6df1SKris Buschelman     nnz     - dim of r, c, and v (number of local nonzero entries of A)
106397b6df1SKris Buschelman     r, c, v - row and col index, matrix values (matrix triples)
107eb9baa12SBarry Smith 
108eb9baa12SBarry Smith   The returned values r, c, and sometimes v are obtained in a single PetscMalloc(). Then in MatDestroy_MUMPS() it is
109eb9baa12SBarry Smith   freed with PetscFree((mumps->irn);  This is not ideal code, the fact that v is ONLY sometimes part of mumps->irn means
110eb9baa12SBarry Smith   that the PetscMalloc() cannot easily be replaced with a PetscMalloc3().
111eb9baa12SBarry Smith 
112397b6df1SKris Buschelman  */
11316ebf90aSShri Abhyankar 
11416ebf90aSShri Abhyankar #undef __FUNCT__
11516ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqaij"
116bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
117b24902e0SBarry Smith {
118185f6596SHong Zhang   const PetscInt *ai,*aj,*ajj,M=A->rmap->n;
11967877ebaSShri Abhyankar   PetscInt       nz,rnz,i,j;
120dfbe8321SBarry Smith   PetscErrorCode ierr;
121c1490034SHong Zhang   PetscInt       *row,*col;
12216ebf90aSShri Abhyankar   Mat_SeqAIJ     *aa=(Mat_SeqAIJ*)A->data;
123397b6df1SKris Buschelman 
124397b6df1SKris Buschelman   PetscFunctionBegin;
12516ebf90aSShri Abhyankar   *v=aa->a;
126bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
1272205254eSKarl Rupp     nz   = aa->nz;
1282205254eSKarl Rupp     ai   = aa->i;
1292205254eSKarl Rupp     aj   = aa->j;
13016ebf90aSShri Abhyankar     *nnz = nz;
131785e854fSJed Brown     ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr);
132185f6596SHong Zhang     col  = row + nz;
133185f6596SHong Zhang 
13416ebf90aSShri Abhyankar     nz = 0;
13516ebf90aSShri Abhyankar     for (i=0; i<M; i++) {
13616ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
13767877ebaSShri Abhyankar       ajj = aj + ai[i];
13867877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
13967877ebaSShri Abhyankar         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
14016ebf90aSShri Abhyankar       }
14116ebf90aSShri Abhyankar     }
14216ebf90aSShri Abhyankar     *r = row; *c = col;
14316ebf90aSShri Abhyankar   }
14416ebf90aSShri Abhyankar   PetscFunctionReturn(0);
14516ebf90aSShri Abhyankar }
146397b6df1SKris Buschelman 
14716ebf90aSShri Abhyankar #undef __FUNCT__
14867877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqbaij_seqaij"
149bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
15067877ebaSShri Abhyankar {
15167877ebaSShri Abhyankar   Mat_SeqBAIJ    *aa=(Mat_SeqBAIJ*)A->data;
15233d57670SJed Brown   const PetscInt *ai,*aj,*ajj,bs2 = aa->bs2;
15333d57670SJed Brown   PetscInt       bs,M,nz,idx=0,rnz,i,j,k,m;
15467877ebaSShri Abhyankar   PetscErrorCode ierr;
15567877ebaSShri Abhyankar   PetscInt       *row,*col;
15667877ebaSShri Abhyankar 
15767877ebaSShri Abhyankar   PetscFunctionBegin;
15833d57670SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
15933d57670SJed Brown   M = A->rmap->N/bs;
160cf3759fdSShri Abhyankar   *v = aa->a;
161bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
162cf3759fdSShri Abhyankar     ai   = aa->i; aj = aa->j;
16367877ebaSShri Abhyankar     nz   = bs2*aa->nz;
16467877ebaSShri Abhyankar     *nnz = nz;
165785e854fSJed Brown     ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr);
166185f6596SHong Zhang     col  = row + nz;
167185f6596SHong Zhang 
16867877ebaSShri Abhyankar     for (i=0; i<M; i++) {
16967877ebaSShri Abhyankar       ajj = aj + ai[i];
17067877ebaSShri Abhyankar       rnz = ai[i+1] - ai[i];
17167877ebaSShri Abhyankar       for (k=0; k<rnz; k++) {
17267877ebaSShri Abhyankar         for (j=0; j<bs; j++) {
17367877ebaSShri Abhyankar           for (m=0; m<bs; m++) {
17467877ebaSShri Abhyankar             row[idx]   = i*bs + m + shift;
175cf3759fdSShri Abhyankar             col[idx++] = bs*(ajj[k]) + j + shift;
17667877ebaSShri Abhyankar           }
17767877ebaSShri Abhyankar         }
17867877ebaSShri Abhyankar       }
17967877ebaSShri Abhyankar     }
180cf3759fdSShri Abhyankar     *r = row; *c = col;
18167877ebaSShri Abhyankar   }
18267877ebaSShri Abhyankar   PetscFunctionReturn(0);
18367877ebaSShri Abhyankar }
18467877ebaSShri Abhyankar 
18567877ebaSShri Abhyankar #undef __FUNCT__
18616ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqsbaij_seqsbaij"
187bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
18816ebf90aSShri Abhyankar {
18967877ebaSShri Abhyankar   const PetscInt *ai, *aj,*ajj,M=A->rmap->n;
19067877ebaSShri Abhyankar   PetscInt       nz,rnz,i,j;
19116ebf90aSShri Abhyankar   PetscErrorCode ierr;
19216ebf90aSShri Abhyankar   PetscInt       *row,*col;
19316ebf90aSShri Abhyankar   Mat_SeqSBAIJ   *aa=(Mat_SeqSBAIJ*)A->data;
19416ebf90aSShri Abhyankar 
19516ebf90aSShri Abhyankar   PetscFunctionBegin;
196882afa5aSHong Zhang   *v = aa->a;
197bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
1982205254eSKarl Rupp     nz   = aa->nz;
1992205254eSKarl Rupp     ai   = aa->i;
2002205254eSKarl Rupp     aj   = aa->j;
2012205254eSKarl Rupp     *v   = aa->a;
20216ebf90aSShri Abhyankar     *nnz = nz;
203785e854fSJed Brown     ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr);
204185f6596SHong Zhang     col  = row + nz;
205185f6596SHong Zhang 
20616ebf90aSShri Abhyankar     nz = 0;
20716ebf90aSShri Abhyankar     for (i=0; i<M; i++) {
20816ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
20967877ebaSShri Abhyankar       ajj = aj + ai[i];
21067877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
21167877ebaSShri Abhyankar         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
21216ebf90aSShri Abhyankar       }
21316ebf90aSShri Abhyankar     }
21416ebf90aSShri Abhyankar     *r = row; *c = col;
21516ebf90aSShri Abhyankar   }
21616ebf90aSShri Abhyankar   PetscFunctionReturn(0);
21716ebf90aSShri Abhyankar }
21816ebf90aSShri Abhyankar 
21916ebf90aSShri Abhyankar #undef __FUNCT__
22016ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqsbaij"
221bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
22216ebf90aSShri Abhyankar {
22367877ebaSShri Abhyankar   const PetscInt    *ai,*aj,*ajj,*adiag,M=A->rmap->n;
22467877ebaSShri Abhyankar   PetscInt          nz,rnz,i,j;
22567877ebaSShri Abhyankar   const PetscScalar *av,*v1;
22616ebf90aSShri Abhyankar   PetscScalar       *val;
22716ebf90aSShri Abhyankar   PetscErrorCode    ierr;
22816ebf90aSShri Abhyankar   PetscInt          *row,*col;
229829b1710SHong Zhang   Mat_SeqAIJ        *aa=(Mat_SeqAIJ*)A->data;
23016ebf90aSShri Abhyankar 
23116ebf90aSShri Abhyankar   PetscFunctionBegin;
23216ebf90aSShri Abhyankar   ai   =aa->i; aj=aa->j;av=aa->a;
23316ebf90aSShri Abhyankar   adiag=aa->diag;
234bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
235829b1710SHong Zhang     /* count nz in the uppper triangular part of A */
236829b1710SHong Zhang     nz = 0;
237829b1710SHong Zhang     for (i=0; i<M; i++) nz += ai[i+1] - adiag[i];
23816ebf90aSShri Abhyankar     *nnz = nz;
239829b1710SHong Zhang 
240185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
241185f6596SHong Zhang     col  = row + nz;
242185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
243185f6596SHong Zhang 
24416ebf90aSShri Abhyankar     nz = 0;
24516ebf90aSShri Abhyankar     for (i=0; i<M; i++) {
24616ebf90aSShri Abhyankar       rnz = ai[i+1] - adiag[i];
24767877ebaSShri Abhyankar       ajj = aj + adiag[i];
248cf3759fdSShri Abhyankar       v1  = av + adiag[i];
24967877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
25067877ebaSShri Abhyankar         row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
25116ebf90aSShri Abhyankar       }
25216ebf90aSShri Abhyankar     }
25316ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
254397b6df1SKris Buschelman   } else {
25516ebf90aSShri Abhyankar     nz = 0; val = *v;
25616ebf90aSShri Abhyankar     for (i=0; i <M; i++) {
25716ebf90aSShri Abhyankar       rnz = ai[i+1] - adiag[i];
25867877ebaSShri Abhyankar       ajj = aj + adiag[i];
25967877ebaSShri Abhyankar       v1  = av + adiag[i];
26067877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
26167877ebaSShri Abhyankar         val[nz++] = v1[j];
26216ebf90aSShri Abhyankar       }
26316ebf90aSShri Abhyankar     }
26416ebf90aSShri Abhyankar   }
26516ebf90aSShri Abhyankar   PetscFunctionReturn(0);
26616ebf90aSShri Abhyankar }
26716ebf90aSShri Abhyankar 
26816ebf90aSShri Abhyankar #undef __FUNCT__
26916ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpisbaij_mpisbaij"
270bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
27116ebf90aSShri Abhyankar {
27216ebf90aSShri Abhyankar   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
27316ebf90aSShri Abhyankar   PetscErrorCode    ierr;
27416ebf90aSShri Abhyankar   PetscInt          rstart,nz,i,j,jj,irow,countA,countB;
27516ebf90aSShri Abhyankar   PetscInt          *row,*col;
27616ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
27716ebf90aSShri Abhyankar   PetscScalar       *val;
278397b6df1SKris Buschelman   Mat_MPISBAIJ      *mat = (Mat_MPISBAIJ*)A->data;
279397b6df1SKris Buschelman   Mat_SeqSBAIJ      *aa  = (Mat_SeqSBAIJ*)(mat->A)->data;
280397b6df1SKris Buschelman   Mat_SeqBAIJ       *bb  = (Mat_SeqBAIJ*)(mat->B)->data;
28116ebf90aSShri Abhyankar 
28216ebf90aSShri Abhyankar   PetscFunctionBegin;
283d0f46423SBarry Smith   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
284397b6df1SKris Buschelman   av=aa->a; bv=bb->a;
285397b6df1SKris Buschelman 
2862205254eSKarl Rupp   garray = mat->garray;
2872205254eSKarl Rupp 
288bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
28916ebf90aSShri Abhyankar     nz   = aa->nz + bb->nz;
29016ebf90aSShri Abhyankar     *nnz = nz;
291185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
292185f6596SHong Zhang     col  = row + nz;
293185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
294185f6596SHong Zhang 
295397b6df1SKris Buschelman     *r = row; *c = col; *v = val;
296397b6df1SKris Buschelman   } else {
297397b6df1SKris Buschelman     row = *r; col = *c; val = *v;
298397b6df1SKris Buschelman   }
299397b6df1SKris Buschelman 
300028e57e8SHong Zhang   jj = 0; irow = rstart;
301397b6df1SKris Buschelman   for (i=0; i<m; i++) {
302397b6df1SKris Buschelman     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
303397b6df1SKris Buschelman     countA = ai[i+1] - ai[i];
304397b6df1SKris Buschelman     countB = bi[i+1] - bi[i];
305397b6df1SKris Buschelman     bjj    = bj + bi[i];
30616ebf90aSShri Abhyankar     v1     = av + ai[i];
30716ebf90aSShri Abhyankar     v2     = bv + bi[i];
308397b6df1SKris Buschelman 
309397b6df1SKris Buschelman     /* A-part */
310397b6df1SKris Buschelman     for (j=0; j<countA; j++) {
311bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
312397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
313397b6df1SKris Buschelman       }
31416ebf90aSShri Abhyankar       val[jj++] = v1[j];
315397b6df1SKris Buschelman     }
31616ebf90aSShri Abhyankar 
31716ebf90aSShri Abhyankar     /* B-part */
31816ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
319bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
320397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
321397b6df1SKris Buschelman       }
32216ebf90aSShri Abhyankar       val[jj++] = v2[j];
32316ebf90aSShri Abhyankar     }
32416ebf90aSShri Abhyankar     irow++;
32516ebf90aSShri Abhyankar   }
32616ebf90aSShri Abhyankar   PetscFunctionReturn(0);
32716ebf90aSShri Abhyankar }
32816ebf90aSShri Abhyankar 
32916ebf90aSShri Abhyankar #undef __FUNCT__
33016ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpiaij"
331bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
33216ebf90aSShri Abhyankar {
33316ebf90aSShri Abhyankar   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
33416ebf90aSShri Abhyankar   PetscErrorCode    ierr;
33516ebf90aSShri Abhyankar   PetscInt          rstart,nz,i,j,jj,irow,countA,countB;
33616ebf90aSShri Abhyankar   PetscInt          *row,*col;
33716ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
33816ebf90aSShri Abhyankar   PetscScalar       *val;
33916ebf90aSShri Abhyankar   Mat_MPIAIJ        *mat = (Mat_MPIAIJ*)A->data;
34016ebf90aSShri Abhyankar   Mat_SeqAIJ        *aa  = (Mat_SeqAIJ*)(mat->A)->data;
34116ebf90aSShri Abhyankar   Mat_SeqAIJ        *bb  = (Mat_SeqAIJ*)(mat->B)->data;
34216ebf90aSShri Abhyankar 
34316ebf90aSShri Abhyankar   PetscFunctionBegin;
34416ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
34516ebf90aSShri Abhyankar   av=aa->a; bv=bb->a;
34616ebf90aSShri Abhyankar 
3472205254eSKarl Rupp   garray = mat->garray;
3482205254eSKarl Rupp 
349bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
35016ebf90aSShri Abhyankar     nz   = aa->nz + bb->nz;
35116ebf90aSShri Abhyankar     *nnz = nz;
352185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
353185f6596SHong Zhang     col  = row + nz;
354185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
355185f6596SHong Zhang 
35616ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
35716ebf90aSShri Abhyankar   } else {
35816ebf90aSShri Abhyankar     row = *r; col = *c; val = *v;
35916ebf90aSShri Abhyankar   }
36016ebf90aSShri Abhyankar 
36116ebf90aSShri Abhyankar   jj = 0; irow = rstart;
36216ebf90aSShri Abhyankar   for (i=0; i<m; i++) {
36316ebf90aSShri Abhyankar     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
36416ebf90aSShri Abhyankar     countA = ai[i+1] - ai[i];
36516ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
36616ebf90aSShri Abhyankar     bjj    = bj + bi[i];
36716ebf90aSShri Abhyankar     v1     = av + ai[i];
36816ebf90aSShri Abhyankar     v2     = bv + bi[i];
36916ebf90aSShri Abhyankar 
37016ebf90aSShri Abhyankar     /* A-part */
37116ebf90aSShri Abhyankar     for (j=0; j<countA; j++) {
372bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
37316ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
37416ebf90aSShri Abhyankar       }
37516ebf90aSShri Abhyankar       val[jj++] = v1[j];
37616ebf90aSShri Abhyankar     }
37716ebf90aSShri Abhyankar 
37816ebf90aSShri Abhyankar     /* B-part */
37916ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
380bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
38116ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
38216ebf90aSShri Abhyankar       }
38316ebf90aSShri Abhyankar       val[jj++] = v2[j];
38416ebf90aSShri Abhyankar     }
38516ebf90aSShri Abhyankar     irow++;
38616ebf90aSShri Abhyankar   }
38716ebf90aSShri Abhyankar   PetscFunctionReturn(0);
38816ebf90aSShri Abhyankar }
38916ebf90aSShri Abhyankar 
39016ebf90aSShri Abhyankar #undef __FUNCT__
39167877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpibaij_mpiaij"
392bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
39367877ebaSShri Abhyankar {
39467877ebaSShri Abhyankar   Mat_MPIBAIJ       *mat    = (Mat_MPIBAIJ*)A->data;
39567877ebaSShri Abhyankar   Mat_SeqBAIJ       *aa     = (Mat_SeqBAIJ*)(mat->A)->data;
39667877ebaSShri Abhyankar   Mat_SeqBAIJ       *bb     = (Mat_SeqBAIJ*)(mat->B)->data;
39767877ebaSShri Abhyankar   const PetscInt    *ai     = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj;
398d985c460SShri Abhyankar   const PetscInt    *garray = mat->garray,mbs=mat->mbs,rstart=A->rmap->rstart;
39933d57670SJed Brown   const PetscInt    bs2=mat->bs2;
40067877ebaSShri Abhyankar   PetscErrorCode    ierr;
40133d57670SJed Brown   PetscInt          bs,nz,i,j,k,n,jj,irow,countA,countB,idx;
40267877ebaSShri Abhyankar   PetscInt          *row,*col;
40367877ebaSShri Abhyankar   const PetscScalar *av=aa->a, *bv=bb->a,*v1,*v2;
40467877ebaSShri Abhyankar   PetscScalar       *val;
40567877ebaSShri Abhyankar 
40667877ebaSShri Abhyankar   PetscFunctionBegin;
40733d57670SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
408bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
40967877ebaSShri Abhyankar     nz   = bs2*(aa->nz + bb->nz);
41067877ebaSShri Abhyankar     *nnz = nz;
411185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
412185f6596SHong Zhang     col  = row + nz;
413185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
414185f6596SHong Zhang 
41567877ebaSShri Abhyankar     *r = row; *c = col; *v = val;
41667877ebaSShri Abhyankar   } else {
41767877ebaSShri Abhyankar     row = *r; col = *c; val = *v;
41867877ebaSShri Abhyankar   }
41967877ebaSShri Abhyankar 
420d985c460SShri Abhyankar   jj = 0; irow = rstart;
42167877ebaSShri Abhyankar   for (i=0; i<mbs; i++) {
42267877ebaSShri Abhyankar     countA = ai[i+1] - ai[i];
42367877ebaSShri Abhyankar     countB = bi[i+1] - bi[i];
42467877ebaSShri Abhyankar     ajj    = aj + ai[i];
42567877ebaSShri Abhyankar     bjj    = bj + bi[i];
42667877ebaSShri Abhyankar     v1     = av + bs2*ai[i];
42767877ebaSShri Abhyankar     v2     = bv + bs2*bi[i];
42867877ebaSShri Abhyankar 
42967877ebaSShri Abhyankar     idx = 0;
43067877ebaSShri Abhyankar     /* A-part */
43167877ebaSShri Abhyankar     for (k=0; k<countA; k++) {
43267877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
43367877ebaSShri Abhyankar         for (n=0; n<bs; n++) {
434bccb9932SShri Abhyankar           if (reuse == MAT_INITIAL_MATRIX) {
435d985c460SShri Abhyankar             row[jj] = irow + n + shift;
436d985c460SShri Abhyankar             col[jj] = rstart + bs*ajj[k] + j + shift;
43767877ebaSShri Abhyankar           }
43867877ebaSShri Abhyankar           val[jj++] = v1[idx++];
43967877ebaSShri Abhyankar         }
44067877ebaSShri Abhyankar       }
44167877ebaSShri Abhyankar     }
44267877ebaSShri Abhyankar 
44367877ebaSShri Abhyankar     idx = 0;
44467877ebaSShri Abhyankar     /* B-part */
44567877ebaSShri Abhyankar     for (k=0; k<countB; k++) {
44667877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
44767877ebaSShri Abhyankar         for (n=0; n<bs; n++) {
448bccb9932SShri Abhyankar           if (reuse == MAT_INITIAL_MATRIX) {
449d985c460SShri Abhyankar             row[jj] = irow + n + shift;
450d985c460SShri Abhyankar             col[jj] = bs*garray[bjj[k]] + j + shift;
45167877ebaSShri Abhyankar           }
452d985c460SShri Abhyankar           val[jj++] = v2[idx++];
45367877ebaSShri Abhyankar         }
45467877ebaSShri Abhyankar       }
45567877ebaSShri Abhyankar     }
456d985c460SShri Abhyankar     irow += bs;
45767877ebaSShri Abhyankar   }
45867877ebaSShri Abhyankar   PetscFunctionReturn(0);
45967877ebaSShri Abhyankar }
46067877ebaSShri Abhyankar 
46167877ebaSShri Abhyankar #undef __FUNCT__
46216ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpisbaij"
463bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
46416ebf90aSShri Abhyankar {
46516ebf90aSShri Abhyankar   const PetscInt    *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
46616ebf90aSShri Abhyankar   PetscErrorCode    ierr;
467e0bace9bSHong Zhang   PetscInt          rstart,nz,nza,nzb,i,j,jj,irow,countA,countB;
46816ebf90aSShri Abhyankar   PetscInt          *row,*col;
46916ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
47016ebf90aSShri Abhyankar   PetscScalar       *val;
47116ebf90aSShri Abhyankar   Mat_MPIAIJ        *mat =  (Mat_MPIAIJ*)A->data;
47216ebf90aSShri Abhyankar   Mat_SeqAIJ        *aa  =(Mat_SeqAIJ*)(mat->A)->data;
47316ebf90aSShri Abhyankar   Mat_SeqAIJ        *bb  =(Mat_SeqAIJ*)(mat->B)->data;
47416ebf90aSShri Abhyankar 
47516ebf90aSShri Abhyankar   PetscFunctionBegin;
47616ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j; adiag=aa->diag;
47716ebf90aSShri Abhyankar   bi=bb->i; bj=bb->j; garray = mat->garray;
47816ebf90aSShri Abhyankar   av=aa->a; bv=bb->a;
4792205254eSKarl Rupp 
48016ebf90aSShri Abhyankar   rstart = A->rmap->rstart;
48116ebf90aSShri Abhyankar 
482bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
483e0bace9bSHong Zhang     nza = 0;    /* num of upper triangular entries in mat->A, including diagonals */
484e0bace9bSHong Zhang     nzb = 0;    /* num of upper triangular entries in mat->B */
48516ebf90aSShri Abhyankar     for (i=0; i<m; i++) {
486e0bace9bSHong Zhang       nza   += (ai[i+1] - adiag[i]);
48716ebf90aSShri Abhyankar       countB = bi[i+1] - bi[i];
48816ebf90aSShri Abhyankar       bjj    = bj + bi[i];
489e0bace9bSHong Zhang       for (j=0; j<countB; j++) {
490e0bace9bSHong Zhang         if (garray[bjj[j]] > rstart) nzb++;
491e0bace9bSHong Zhang       }
492e0bace9bSHong Zhang     }
49316ebf90aSShri Abhyankar 
494e0bace9bSHong Zhang     nz   = nza + nzb; /* total nz of upper triangular part of mat */
49516ebf90aSShri Abhyankar     *nnz = nz;
496185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
497185f6596SHong Zhang     col  = row + nz;
498185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
499185f6596SHong Zhang 
50016ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
50116ebf90aSShri Abhyankar   } else {
50216ebf90aSShri Abhyankar     row = *r; col = *c; val = *v;
50316ebf90aSShri Abhyankar   }
50416ebf90aSShri Abhyankar 
50516ebf90aSShri Abhyankar   jj = 0; irow = rstart;
50616ebf90aSShri Abhyankar   for (i=0; i<m; i++) {
50716ebf90aSShri Abhyankar     ajj    = aj + adiag[i];                 /* ptr to the beginning of the diagonal of this row */
50816ebf90aSShri Abhyankar     v1     = av + adiag[i];
50916ebf90aSShri Abhyankar     countA = ai[i+1] - adiag[i];
51016ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
51116ebf90aSShri Abhyankar     bjj    = bj + bi[i];
51216ebf90aSShri Abhyankar     v2     = bv + bi[i];
51316ebf90aSShri Abhyankar 
51416ebf90aSShri Abhyankar     /* A-part */
51516ebf90aSShri Abhyankar     for (j=0; j<countA; j++) {
516bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
51716ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
51816ebf90aSShri Abhyankar       }
51916ebf90aSShri Abhyankar       val[jj++] = v1[j];
52016ebf90aSShri Abhyankar     }
52116ebf90aSShri Abhyankar 
52216ebf90aSShri Abhyankar     /* B-part */
52316ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
52416ebf90aSShri Abhyankar       if (garray[bjj[j]] > rstart) {
525bccb9932SShri Abhyankar         if (reuse == MAT_INITIAL_MATRIX) {
52616ebf90aSShri Abhyankar           row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
52716ebf90aSShri Abhyankar         }
52816ebf90aSShri Abhyankar         val[jj++] = v2[j];
52916ebf90aSShri Abhyankar       }
530397b6df1SKris Buschelman     }
531397b6df1SKris Buschelman     irow++;
532397b6df1SKris Buschelman   }
533397b6df1SKris Buschelman   PetscFunctionReturn(0);
534397b6df1SKris Buschelman }
535397b6df1SKris Buschelman 
536397b6df1SKris Buschelman #undef __FUNCT__
53720be8e61SHong Zhang #define __FUNCT__ "MatGetDiagonal_MUMPS"
53820be8e61SHong Zhang PetscErrorCode MatGetDiagonal_MUMPS(Mat A,Vec v)
53920be8e61SHong Zhang {
54020be8e61SHong Zhang   PetscFunctionBegin;
54120be8e61SHong Zhang   SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type: MUMPS factor");
54220be8e61SHong Zhang   PetscFunctionReturn(0);
54320be8e61SHong Zhang }
54420be8e61SHong Zhang 
54520be8e61SHong Zhang #undef __FUNCT__
5463924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS"
547dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A)
548dfbe8321SBarry Smith {
549a5e57a09SHong Zhang   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
550dfbe8321SBarry Smith   PetscErrorCode ierr;
551b24902e0SBarry Smith 
552397b6df1SKris Buschelman   PetscFunctionBegin;
553a5e57a09SHong Zhang   if (mumps->CleanUpMUMPS) {
554397b6df1SKris Buschelman     /* Terminate instance, deallocate memories */
555a5e57a09SHong Zhang     ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr);
556a5e57a09SHong Zhang     ierr = VecScatterDestroy(&mumps->scat_rhs);CHKERRQ(ierr);
557a5e57a09SHong Zhang     ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
558801fbe65SHong Zhang     ierr = VecDestroy(&mumps->b_seq);CHKERRQ(ierr);
559a5e57a09SHong Zhang     ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr);
560a5e57a09SHong Zhang     ierr = PetscFree(mumps->id.perm_in);CHKERRQ(ierr);
561a5e57a09SHong Zhang     ierr = PetscFree(mumps->irn);CHKERRQ(ierr);
562b34f08ffSHong Zhang     ierr = PetscFree(mumps->info);CHKERRQ(ierr);
5632205254eSKarl Rupp 
564a5e57a09SHong Zhang     mumps->id.job = JOB_END;
565a5e57a09SHong Zhang     PetscMUMPS_c(&mumps->id);
566a5e57a09SHong Zhang     ierr = MPI_Comm_free(&(mumps->comm_mumps));CHKERRQ(ierr);
567397b6df1SKris Buschelman   }
568a5e57a09SHong Zhang   if (mumps->Destroy) {
569a5e57a09SHong Zhang     ierr = (mumps->Destroy)(A);CHKERRQ(ierr);
570bf0cc555SLisandro Dalcin   }
571bf0cc555SLisandro Dalcin   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
572bf0cc555SLisandro Dalcin 
57397969023SHong Zhang   /* clear composed functions */
574bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr);
575bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetIcntl_C",NULL);CHKERRQ(ierr);
576bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetIcntl_C",NULL);CHKERRQ(ierr);
577bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetCntl_C",NULL);CHKERRQ(ierr);
578bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetCntl_C",NULL);CHKERRQ(ierr);
579bc6112feSHong Zhang 
580ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfo_C",NULL);CHKERRQ(ierr);
581ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfog_C",NULL);CHKERRQ(ierr);
582ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfo_C",NULL);CHKERRQ(ierr);
583ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfog_C",NULL);CHKERRQ(ierr);
584397b6df1SKris Buschelman   PetscFunctionReturn(0);
585397b6df1SKris Buschelman }
586397b6df1SKris Buschelman 
587397b6df1SKris Buschelman #undef __FUNCT__
588f6c57405SHong Zhang #define __FUNCT__ "MatSolve_MUMPS"
589b24902e0SBarry Smith PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x)
590b24902e0SBarry Smith {
591a5e57a09SHong Zhang   Mat_MUMPS        *mumps=(Mat_MUMPS*)A->spptr;
592d54de34fSKris Buschelman   PetscScalar      *array;
59367877ebaSShri Abhyankar   Vec              b_seq;
594329ec9b3SHong Zhang   IS               is_iden,is_petsc;
595dfbe8321SBarry Smith   PetscErrorCode   ierr;
596329ec9b3SHong Zhang   PetscInt         i;
597883f2eb9SBarry Smith   static PetscBool cite1 = PETSC_FALSE,cite2 = PETSC_FALSE;
598397b6df1SKris Buschelman 
599397b6df1SKris Buschelman   PetscFunctionBegin;
600883f2eb9SBarry Smith   ierr = PetscCitationsRegister("@article{MUMPS01,\n  author = {P.~R. Amestoy and I.~S. Duff and J.-Y. L'Excellent and J. Koster},\n  title = {A fully asynchronous multifrontal solver using distributed dynamic scheduling},\n  journal = {SIAM Journal on Matrix Analysis and Applications},\n  volume = {23},\n  number = {1},\n  pages = {15--41},\n  year = {2001}\n}\n",&cite1);CHKERRQ(ierr);
601883f2eb9SBarry Smith   ierr = PetscCitationsRegister("@article{MUMPS02,\n  author = {P.~R. Amestoy and A. Guermouche and J.-Y. L'Excellent and S. Pralet},\n  title = {Hybrid scheduling for the parallel solution of linear systems},\n  journal = {Parallel Computing},\n  volume = {32},\n  number = {2},\n  pages = {136--156},\n  year = {2006}\n}\n",&cite2);CHKERRQ(ierr);
602a5e57a09SHong Zhang   mumps->id.nrhs = 1;
603a5e57a09SHong Zhang   b_seq          = mumps->b_seq;
604a5e57a09SHong Zhang   if (mumps->size > 1) {
605329ec9b3SHong Zhang     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
606a5e57a09SHong Zhang     ierr = VecScatterBegin(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
607a5e57a09SHong Zhang     ierr = VecScatterEnd(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
608a5e57a09SHong Zhang     if (!mumps->myid) {ierr = VecGetArray(b_seq,&array);CHKERRQ(ierr);}
609397b6df1SKris Buschelman   } else {  /* size == 1 */
610397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
611397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
612397b6df1SKris Buschelman   }
613a5e57a09SHong Zhang   if (!mumps->myid) { /* define rhs on the host */
614a5e57a09SHong Zhang     mumps->id.nrhs = 1;
615*940cd9d6SSatish Balay     mumps->id.rhs = (MumpsScalar*)array;
616397b6df1SKris Buschelman   }
617397b6df1SKris Buschelman 
618397b6df1SKris Buschelman   /* solve phase */
619329ec9b3SHong Zhang   /*-------------*/
620a5e57a09SHong Zhang   mumps->id.job = JOB_SOLVE;
621a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
622a5e57a09SHong 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));
623397b6df1SKris Buschelman 
624a5e57a09SHong Zhang   if (mumps->size > 1) { /* convert mumps distributed solution to petsc mpi x */
625a5e57a09SHong Zhang     if (mumps->scat_sol && mumps->ICNTL9_pre != mumps->id.ICNTL(9)) {
626a5e57a09SHong Zhang       /* when id.ICNTL(9) changes, the contents of lsol_loc may change (not its size, lsol_loc), recreates scat_sol */
627a5e57a09SHong Zhang       ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
628397b6df1SKris Buschelman     }
629a5e57a09SHong Zhang     if (!mumps->scat_sol) { /* create scatter scat_sol */
630a5e57a09SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,mumps->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */
631a5e57a09SHong Zhang       for (i=0; i<mumps->id.lsol_loc; i++) {
632a5e57a09SHong Zhang         mumps->id.isol_loc[i] -= 1; /* change Fortran style to C style */
633a5e57a09SHong Zhang       }
634a5e57a09SHong Zhang       ierr = ISCreateGeneral(PETSC_COMM_SELF,mumps->id.lsol_loc,mumps->id.isol_loc,PETSC_COPY_VALUES,&is_petsc);CHKERRQ(ierr);  /* to */
635a5e57a09SHong Zhang       ierr = VecScatterCreate(mumps->x_seq,is_iden,x,is_petsc,&mumps->scat_sol);CHKERRQ(ierr);
6366bf464f9SBarry Smith       ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
6376bf464f9SBarry Smith       ierr = ISDestroy(&is_petsc);CHKERRQ(ierr);
6382205254eSKarl Rupp 
639a5e57a09SHong Zhang       mumps->ICNTL9_pre = mumps->id.ICNTL(9); /* save current value of id.ICNTL(9) */
640397b6df1SKris Buschelman     }
641a5e57a09SHong Zhang 
642a5e57a09SHong Zhang     ierr = VecScatterBegin(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
643a5e57a09SHong Zhang     ierr = VecScatterEnd(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
644329ec9b3SHong Zhang   }
645397b6df1SKris Buschelman   PetscFunctionReturn(0);
646397b6df1SKris Buschelman }
647397b6df1SKris Buschelman 
64851d5961aSHong Zhang #undef __FUNCT__
64951d5961aSHong Zhang #define __FUNCT__ "MatSolveTranspose_MUMPS"
65051d5961aSHong Zhang PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x)
65151d5961aSHong Zhang {
652a5e57a09SHong Zhang   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
65351d5961aSHong Zhang   PetscErrorCode ierr;
65451d5961aSHong Zhang 
65551d5961aSHong Zhang   PetscFunctionBegin;
656a5e57a09SHong Zhang   mumps->id.ICNTL(9) = 0;
6570ad0caddSJed Brown   ierr = MatSolve_MUMPS(A,b,x);CHKERRQ(ierr);
658a5e57a09SHong Zhang   mumps->id.ICNTL(9) = 1;
65951d5961aSHong Zhang   PetscFunctionReturn(0);
66051d5961aSHong Zhang }
66151d5961aSHong Zhang 
662e0b74bf9SHong Zhang #undef __FUNCT__
663e0b74bf9SHong Zhang #define __FUNCT__ "MatMatSolve_MUMPS"
664e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X)
665e0b74bf9SHong Zhang {
666bda8bf91SBarry Smith   PetscErrorCode ierr;
667bda8bf91SBarry Smith   PetscBool      flg;
6684e34a73bSHong Zhang   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
669334c5f61SHong Zhang   PetscInt       i,nrhs,M;
6702cd7d884SHong Zhang   PetscScalar    *array,*bray;
671bda8bf91SBarry Smith 
672e0b74bf9SHong Zhang   PetscFunctionBegin;
6730298fd71SBarry Smith   ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
674801fbe65SHong Zhang   if (!flg) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
6750298fd71SBarry Smith   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
676801fbe65SHong Zhang   if (!flg) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
677801fbe65SHong Zhang   if (B->rmap->n != X->rmap->n) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B and X must have same row distribution");
6784e34a73bSHong Zhang 
6792cd7d884SHong Zhang   ierr = MatGetSize(B,&M,&nrhs);CHKERRQ(ierr);
680334c5f61SHong Zhang   mumps->id.nrhs = nrhs;
681334c5f61SHong Zhang   mumps->id.lrhs = M;
6824e34a73bSHong Zhang 
6832cd7d884SHong Zhang   if (mumps->size == 1) {
6842cd7d884SHong Zhang     /* copy B to X */
6852cd7d884SHong Zhang     ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr);
6862cd7d884SHong Zhang     ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr);
6872cd7d884SHong Zhang     for (i=0; i<M*nrhs; i++) array[i] = bray[i];
6882cd7d884SHong Zhang     ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr);
689*940cd9d6SSatish Balay     mumps->id.rhs = (MumpsScalar*)array;
690801fbe65SHong Zhang 
6912cd7d884SHong Zhang     /* solve phase */
6922cd7d884SHong Zhang     /*-------------*/
6932cd7d884SHong Zhang     mumps->id.job = JOB_SOLVE;
6942cd7d884SHong Zhang     PetscMUMPS_c(&mumps->id);
6952cd7d884SHong 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));
6962cd7d884SHong Zhang     ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr);
697334c5f61SHong Zhang   } else {  /*--------- parallel case --------*/
69871aed81dSHong Zhang     PetscInt       lsol_loc,nlsol_loc,*isol_loc,*idx,*iidx,*idxx,*isol_loc_save;
69971aed81dSHong Zhang     PetscScalar    *sol_loc,*sol_loc_save;
700801fbe65SHong Zhang     IS             is_to,is_from;
701334c5f61SHong Zhang     PetscInt       k,proc,j,m;
702801fbe65SHong Zhang     const PetscInt *rstart;
703334c5f61SHong Zhang     Vec            v_mpi,b_seq,x_seq;
704334c5f61SHong Zhang     VecScatter     scat_rhs,scat_sol;
705801fbe65SHong Zhang 
706801fbe65SHong Zhang     /* create x_seq to hold local solution */
70771aed81dSHong Zhang     isol_loc_save = mumps->id.isol_loc; /* save it for MatSovle() */
70871aed81dSHong Zhang     sol_loc_save  = mumps->id.sol_loc;
709801fbe65SHong Zhang 
71071aed81dSHong Zhang     lsol_loc  = mumps->id.INFO(23);
71171aed81dSHong Zhang     nlsol_loc = nrhs*lsol_loc;     /* length of sol_loc */
71271aed81dSHong Zhang     ierr = PetscMalloc2(nlsol_loc,&sol_loc,nlsol_loc,&isol_loc);CHKERRQ(ierr);
713*940cd9d6SSatish Balay     mumps->id.sol_loc = (MumpsScalar*)sol_loc;
714801fbe65SHong Zhang     mumps->id.isol_loc = isol_loc;
715801fbe65SHong Zhang 
71671aed81dSHong Zhang     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,nlsol_loc,sol_loc,&x_seq);CHKERRQ(ierr);
7172cd7d884SHong Zhang 
71874f0fcc7SHong Zhang     /* copy rhs matrix B into vector v_mpi */
719334c5f61SHong Zhang     ierr = MatGetLocalSize(B,&m,NULL);CHKERRQ(ierr);
720801fbe65SHong Zhang     ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr);
72174f0fcc7SHong Zhang     ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)B),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi);CHKERRQ(ierr);
722801fbe65SHong Zhang     ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr);
723801fbe65SHong Zhang 
724334c5f61SHong Zhang     /* scatter v_mpi to b_seq because MUMPS only supports centralized rhs */
72574f0fcc7SHong Zhang     /* idx: maps from k-th index of v_mpi to (i,j)-th global entry of B;
726801fbe65SHong Zhang       iidx: inverse of idx, will be used by scattering xx_seq -> X       */
727801fbe65SHong Zhang     ierr = PetscMalloc2(nrhs*M,&idx,nrhs*M,&iidx);CHKERRQ(ierr);
728801fbe65SHong Zhang     ierr = MatGetOwnershipRanges(B,&rstart);CHKERRQ(ierr);
729801fbe65SHong Zhang     k = 0;
730801fbe65SHong Zhang     for (proc=0; proc<mumps->size; proc++){
731801fbe65SHong Zhang       for (j=0; j<nrhs; j++){
732801fbe65SHong Zhang         for (i=rstart[proc]; i<rstart[proc+1]; i++){
733801fbe65SHong Zhang           iidx[j*M + i] = k;
734801fbe65SHong Zhang           idx[k++]      = j*M + i;
735801fbe65SHong Zhang         }
736801fbe65SHong Zhang       }
7372cd7d884SHong Zhang     }
7382cd7d884SHong Zhang 
739801fbe65SHong Zhang     if (!mumps->myid) {
740334c5f61SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,nrhs*M,&b_seq);CHKERRQ(ierr);
741801fbe65SHong Zhang       ierr = ISCreateGeneral(PETSC_COMM_SELF,nrhs*M,idx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr);
742801fbe65SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,nrhs*M,0,1,&is_from);CHKERRQ(ierr);
743801fbe65SHong Zhang     } else {
744334c5f61SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&b_seq);CHKERRQ(ierr);
745801fbe65SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_to);CHKERRQ(ierr);
746801fbe65SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_from);CHKERRQ(ierr);
747801fbe65SHong Zhang     }
748334c5f61SHong Zhang     ierr = VecScatterCreate(v_mpi,is_from,b_seq,is_to,&scat_rhs);CHKERRQ(ierr);
749334c5f61SHong Zhang     ierr = VecScatterBegin(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
750801fbe65SHong Zhang     ierr = ISDestroy(&is_to);CHKERRQ(ierr);
751801fbe65SHong Zhang     ierr = ISDestroy(&is_from);CHKERRQ(ierr);
752334c5f61SHong Zhang     ierr = VecScatterEnd(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
753801fbe65SHong Zhang 
754801fbe65SHong Zhang     if (!mumps->myid) { /* define rhs on the host */
755334c5f61SHong Zhang       ierr = VecGetArray(b_seq,&bray);CHKERRQ(ierr);
756*940cd9d6SSatish Balay       mumps->id.rhs = (MumpsScalar*)bray;
757334c5f61SHong Zhang       ierr = VecRestoreArray(b_seq,&bray);CHKERRQ(ierr);
758801fbe65SHong Zhang     }
759801fbe65SHong Zhang 
760801fbe65SHong Zhang     /* solve phase */
761801fbe65SHong Zhang     /*-------------*/
762801fbe65SHong Zhang     mumps->id.job = JOB_SOLVE;
763801fbe65SHong Zhang     PetscMUMPS_c(&mumps->id);
764801fbe65SHong 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));
765801fbe65SHong Zhang 
766334c5f61SHong Zhang     /* scatter mumps distributed solution to petsc vector v_mpi, which shares local arrays with solution matrix X */
76774f0fcc7SHong Zhang     ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr);
76874f0fcc7SHong Zhang     ierr = VecPlaceArray(v_mpi,array);CHKERRQ(ierr);
769801fbe65SHong Zhang 
770334c5f61SHong Zhang     /* create scatter scat_sol */
77171aed81dSHong Zhang     ierr = PetscMalloc1(nlsol_loc,&idxx);CHKERRQ(ierr);
77271aed81dSHong Zhang     ierr = ISCreateStride(PETSC_COMM_SELF,nlsol_loc,0,1,&is_from);CHKERRQ(ierr);
77371aed81dSHong Zhang     for (i=0; i<lsol_loc; i++) {
774334c5f61SHong Zhang       isol_loc[i] -= 1; /* change Fortran style to C style */
775334c5f61SHong Zhang       idxx[i] = iidx[isol_loc[i]];
776801fbe65SHong Zhang       for (j=1; j<nrhs; j++){
777334c5f61SHong Zhang         idxx[j*lsol_loc+i] = iidx[isol_loc[i]+j*M];
778801fbe65SHong Zhang       }
779801fbe65SHong Zhang     }
78071aed81dSHong Zhang     ierr = ISCreateGeneral(PETSC_COMM_SELF,nlsol_loc,idxx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr);
781334c5f61SHong Zhang     ierr = VecScatterCreate(x_seq,is_from,v_mpi,is_to,&scat_sol);CHKERRQ(ierr);
782334c5f61SHong Zhang     ierr = VecScatterBegin(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
783801fbe65SHong Zhang     ierr = ISDestroy(&is_from);CHKERRQ(ierr);
784801fbe65SHong Zhang     ierr = ISDestroy(&is_to);CHKERRQ(ierr);
785334c5f61SHong Zhang     ierr = VecScatterEnd(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
786801fbe65SHong Zhang     ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr);
78771aed81dSHong Zhang 
78871aed81dSHong Zhang     /* free spaces */
78971aed81dSHong Zhang     mumps->id.sol_loc = sol_loc_save;
79071aed81dSHong Zhang     mumps->id.isol_loc = isol_loc_save;
79171aed81dSHong Zhang 
79271aed81dSHong Zhang     ierr = PetscFree2(sol_loc,isol_loc);CHKERRQ(ierr);
793801fbe65SHong Zhang     ierr = PetscFree2(idx,iidx);CHKERRQ(ierr);
794801fbe65SHong Zhang     ierr = PetscFree(idxx);CHKERRQ(ierr);
79571aed81dSHong Zhang     ierr = VecDestroy(&x_seq);CHKERRQ(ierr);
79674f0fcc7SHong Zhang     ierr = VecDestroy(&v_mpi);CHKERRQ(ierr);
797334c5f61SHong Zhang     ierr = VecDestroy(&b_seq);CHKERRQ(ierr);
798334c5f61SHong Zhang     ierr = VecScatterDestroy(&scat_rhs);CHKERRQ(ierr);
799334c5f61SHong Zhang     ierr = VecScatterDestroy(&scat_sol);CHKERRQ(ierr);
800801fbe65SHong Zhang   }
801e0b74bf9SHong Zhang   PetscFunctionReturn(0);
802e0b74bf9SHong Zhang }
803e0b74bf9SHong Zhang 
804ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX)
805a58c3f20SHong Zhang /*
806a58c3f20SHong Zhang   input:
807a58c3f20SHong Zhang    F:        numeric factor
808a58c3f20SHong Zhang   output:
809a58c3f20SHong Zhang    nneg:     total number of negative pivots
810a58c3f20SHong Zhang    nzero:    0
811a58c3f20SHong Zhang    npos:     (global dimension of F) - nneg
812a58c3f20SHong Zhang */
813a58c3f20SHong Zhang 
814a58c3f20SHong Zhang #undef __FUNCT__
815a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
816dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
817a58c3f20SHong Zhang {
818a5e57a09SHong Zhang   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
819dfbe8321SBarry Smith   PetscErrorCode ierr;
820c1490034SHong Zhang   PetscMPIInt    size;
821a58c3f20SHong Zhang 
822a58c3f20SHong Zhang   PetscFunctionBegin;
823ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)F),&size);CHKERRQ(ierr);
824bcb30aebSHong 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 */
825a5e57a09SHong 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));
826ed85ac9fSHong Zhang 
827710ac8efSHong Zhang   if (nneg) *nneg = mumps->id.INFOG(12);
828ed85ac9fSHong Zhang   if (nzero || npos) {
829ed85ac9fSHong Zhang     if (mumps->id.ICNTL(24) != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"-mat_mumps_icntl_24 must be set as 1 for null pivot row detection");
830710ac8efSHong Zhang     if (nzero) *nzero = mumps->id.INFOG(28);
831710ac8efSHong Zhang     if (npos) *npos   = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28));
832a58c3f20SHong Zhang   }
833a58c3f20SHong Zhang   PetscFunctionReturn(0);
834a58c3f20SHong Zhang }
835ace3df97SHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */
836a58c3f20SHong Zhang 
837397b6df1SKris Buschelman #undef __FUNCT__
838f6c57405SHong Zhang #define __FUNCT__ "MatFactorNumeric_MUMPS"
8390481f469SBarry Smith PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
840af281ebdSHong Zhang {
841a5e57a09SHong Zhang   Mat_MUMPS      *mumps =(Mat_MUMPS*)(F)->spptr;
8426849ba73SBarry Smith   PetscErrorCode ierr;
843e09efc27SHong Zhang   Mat            F_diag;
844ace3abfcSBarry Smith   PetscBool      isMPIAIJ;
845397b6df1SKris Buschelman 
846397b6df1SKris Buschelman   PetscFunctionBegin;
847a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
848397b6df1SKris Buschelman 
849397b6df1SKris Buschelman   /* numerical factorization phase */
850329ec9b3SHong Zhang   /*-------------------------------*/
851a5e57a09SHong Zhang   mumps->id.job = JOB_FACTNUMERIC;
8524e34a73bSHong Zhang   if (!mumps->id.ICNTL(18)) { /* A is centralized */
853a5e57a09SHong Zhang     if (!mumps->myid) {
854*940cd9d6SSatish Balay       mumps->id.a = (MumpsScalar*)mumps->val;
855397b6df1SKris Buschelman     }
856397b6df1SKris Buschelman   } else {
857*940cd9d6SSatish Balay     mumps->id.a_loc = (MumpsScalar*)mumps->val;
858397b6df1SKris Buschelman   }
859a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
860a5e57a09SHong Zhang   if (mumps->id.INFOG(1) < 0) {
861151787a6SHong Zhang     if (mumps->id.INFO(1) == -13) {
862151787a6SHong Zhang       if (mumps->id.INFO(2) < 0) {
863151787a6SHong Zhang         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",-mumps->id.INFO(2));
864151787a6SHong Zhang       } else {
865151787a6SHong Zhang         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d bytes\n",mumps->id.INFO(2));
866151787a6SHong Zhang       }
867151787a6SHong 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));
868397b6df1SKris Buschelman   }
869a5e57a09SHong 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));
870397b6df1SKris Buschelman 
871dcd589f8SShri Abhyankar   (F)->assembled      = PETSC_TRUE;
872a5e57a09SHong Zhang   mumps->matstruc     = SAME_NONZERO_PATTERN;
873a5e57a09SHong Zhang   mumps->CleanUpMUMPS = PETSC_TRUE;
87467877ebaSShri Abhyankar 
875a5e57a09SHong Zhang   if (mumps->size > 1) {
87667877ebaSShri Abhyankar     PetscInt    lsol_loc;
87767877ebaSShri Abhyankar     PetscScalar *sol_loc;
8782205254eSKarl Rupp 
879c2093ab7SHong Zhang     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr);
880c2093ab7SHong Zhang     if (isMPIAIJ) F_diag = ((Mat_MPIAIJ*)(F)->data)->A;
881c2093ab7SHong Zhang     else F_diag = ((Mat_MPISBAIJ*)(F)->data)->A;
882c2093ab7SHong Zhang     F_diag->assembled = PETSC_TRUE;
883c2093ab7SHong Zhang 
884c2093ab7SHong Zhang     /* distributed solution; Create x_seq=sol_loc for repeated use */
885c2093ab7SHong Zhang     if (mumps->x_seq) {
886c2093ab7SHong Zhang       ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
887c2093ab7SHong Zhang       ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr);
888c2093ab7SHong Zhang       ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr);
889c2093ab7SHong Zhang     }
890a5e57a09SHong Zhang     lsol_loc = mumps->id.INFO(23); /* length of sol_loc */
891dcca6d9dSJed Brown     ierr = PetscMalloc2(lsol_loc,&sol_loc,lsol_loc,&mumps->id.isol_loc);CHKERRQ(ierr);
892a5e57a09SHong Zhang     mumps->id.lsol_loc = lsol_loc;
893*940cd9d6SSatish Balay     mumps->id.sol_loc = (MumpsScalar*)sol_loc;
894a5e57a09SHong Zhang     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq);CHKERRQ(ierr);
89567877ebaSShri Abhyankar   }
896397b6df1SKris Buschelman   PetscFunctionReturn(0);
897397b6df1SKris Buschelman }
898397b6df1SKris Buschelman 
8999a2535b5SHong Zhang /* Sets MUMPS options from the options database */
900dcd589f8SShri Abhyankar #undef __FUNCT__
9019a2535b5SHong Zhang #define __FUNCT__ "PetscSetMUMPSFromOptions"
9029a2535b5SHong Zhang PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A)
903dcd589f8SShri Abhyankar {
9049a2535b5SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
905dcd589f8SShri Abhyankar   PetscErrorCode ierr;
906b34f08ffSHong Zhang   PetscInt       icntl,info[40],i,ninfo=40;
907ace3abfcSBarry Smith   PetscBool      flg;
908dcd589f8SShri Abhyankar 
909dcd589f8SShri Abhyankar   PetscFunctionBegin;
910ce94432eSBarry Smith   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
9119a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);CHKERRQ(ierr);
9129a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(1) = icntl;
9139a2535b5SHong 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);
9149a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(2) = icntl;
9159a2535b5SHong 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);
9169a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(3) = icntl;
917dcd589f8SShri Abhyankar 
9189a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",mumps->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
9199a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(4) = icntl;
9209a2535b5SHong Zhang   if (mumps->id.ICNTL(4) || PetscLogPrintInfo) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */
9219a2535b5SHong Zhang 
922d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_6","ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7)","None",mumps->id.ICNTL(6),&icntl,&flg);CHKERRQ(ierr);
9239a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(6) = icntl;
9249a2535b5SHong Zhang 
925d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis","None",mumps->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
926dcd589f8SShri Abhyankar   if (flg) {
9272205254eSKarl 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");
9282205254eSKarl Rupp     else mumps->id.ICNTL(7) = icntl;
929dcd589f8SShri Abhyankar   }
930e0b74bf9SHong Zhang 
9310298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_8","ICNTL(8): scaling strategy (-2 to 8 or 77)","None",mumps->id.ICNTL(8),&mumps->id.ICNTL(8),NULL);CHKERRQ(ierr);
932d341cd04SHong Zhang   /* ierr = PetscOptionsInt("-mat_mumps_icntl_9","ICNTL(9): computes the solution using A or A^T","None",mumps->id.ICNTL(9),&mumps->id.ICNTL(9),NULL);CHKERRQ(ierr); handled by MatSolveTranspose_MUMPS() */
9330298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",mumps->id.ICNTL(10),&mumps->id.ICNTL(10),NULL);CHKERRQ(ierr);
934d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_11","ICNTL(11): statistics related to an error analysis (via -ksp_view)","None",mumps->id.ICNTL(11),&mumps->id.ICNTL(11),NULL);CHKERRQ(ierr);
935d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): an ordering strategy for symmetric matrices (0 to 3)","None",mumps->id.ICNTL(12),&mumps->id.ICNTL(12),NULL);CHKERRQ(ierr);
936d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting","None",mumps->id.ICNTL(13),&mumps->id.ICNTL(13),NULL);CHKERRQ(ierr);
937d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage increase in the estimated working space","None",mumps->id.ICNTL(14),&mumps->id.ICNTL(14),NULL);CHKERRQ(ierr);
938d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): computes the Schur complement","None",mumps->id.ICNTL(19),&mumps->id.ICNTL(19),NULL);CHKERRQ(ierr);
9394e34a73bSHong Zhang   /* ierr = PetscOptionsInt("-mat_mumps_icntl_20","ICNTL(20): the format (dense or sparse) of the right-hand sides","None",mumps->id.ICNTL(20),&mumps->id.ICNTL(20),NULL);CHKERRQ(ierr); -- sparse rhs is not supported in PETSc API */
940d341cd04SHong Zhang   /* ierr = PetscOptionsInt("-mat_mumps_icntl_21","ICNTL(21): the distribution (centralized or distributed) of the solution vectors","None",mumps->id.ICNTL(21),&mumps->id.ICNTL(21),NULL);CHKERRQ(ierr); we only use distributed solution vector */
9419a2535b5SHong Zhang 
942d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_22","ICNTL(22): in-core/out-of-core factorization and solve (0 or 1)","None",mumps->id.ICNTL(22),&mumps->id.ICNTL(22),NULL);CHKERRQ(ierr);
9430298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_23","ICNTL(23): max size of the working memory (MB) that can allocate per processor","None",mumps->id.ICNTL(23),&mumps->id.ICNTL(23),NULL);CHKERRQ(ierr);
9440298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_24","ICNTL(24): detection of null pivot rows (0 or 1)","None",mumps->id.ICNTL(24),&mumps->id.ICNTL(24),NULL);CHKERRQ(ierr);
9459a2535b5SHong Zhang   if (mumps->id.ICNTL(24)) {
9469a2535b5SHong Zhang     mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */
947d7ebd59bSHong Zhang   }
948d7ebd59bSHong Zhang 
949d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_25","ICNTL(25): compute a solution of a deficient matrix and a null space basis","None",mumps->id.ICNTL(25),&mumps->id.ICNTL(25),NULL);CHKERRQ(ierr);
950d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_26","ICNTL(26): drives the solution phase if a Schur complement matrix","None",mumps->id.ICNTL(26),&mumps->id.ICNTL(26),NULL);CHKERRQ(ierr);
9512cd7d884SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_27","ICNTL(27): the blocking size for multiple right-hand sides","None",mumps->id.ICNTL(27),&mumps->id.ICNTL(27),NULL);CHKERRQ(ierr);
9520298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_28","ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering","None",mumps->id.ICNTL(28),&mumps->id.ICNTL(28),NULL);CHKERRQ(ierr);
953d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_29","ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis","None",mumps->id.ICNTL(29),&mumps->id.ICNTL(29),NULL);CHKERRQ(ierr);
9540298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_30","ICNTL(30): compute user-specified set of entries in inv(A)","None",mumps->id.ICNTL(30),&mumps->id.ICNTL(30),NULL);CHKERRQ(ierr);
955d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_31","ICNTL(31): indicates which factors may be discarded during factorization","None",mumps->id.ICNTL(31),&mumps->id.ICNTL(31),NULL);CHKERRQ(ierr);
9564e34a73bSHong Zhang   /* ierr = PetscOptionsInt("-mat_mumps_icntl_32","ICNTL(32): performs the forward elemination of the right-hand sides during factorization","None",mumps->id.ICNTL(32),&mumps->id.ICNTL(32),NULL);CHKERRQ(ierr);  -- not supported by PETSc API */
9570298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),NULL);CHKERRQ(ierr);
958dcd589f8SShri Abhyankar 
9590298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),NULL);CHKERRQ(ierr);
9600298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",mumps->id.CNTL(2),&mumps->id.CNTL(2),NULL);CHKERRQ(ierr);
9610298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),NULL);CHKERRQ(ierr);
9620298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",mumps->id.CNTL(4),&mumps->id.CNTL(4),NULL);CHKERRQ(ierr);
9630298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",mumps->id.CNTL(5),&mumps->id.CNTL(5),NULL);CHKERRQ(ierr);
964e5bb22a1SHong Zhang 
9650298fd71SBarry Smith   ierr = PetscOptionsString("-mat_mumps_ooc_tmpdir", "out of core directory", "None", mumps->id.ooc_tmpdir, mumps->id.ooc_tmpdir, 256, NULL);
966b34f08ffSHong Zhang 
967b34f08ffSHong Zhang   ierr = PetscOptionsGetIntArray(NULL,"-mat_mumps_view_info",info,&ninfo,NULL);CHKERRQ(ierr); /* why does not show with '-help'??? */
968b34f08ffSHong Zhang   if (ninfo) {
969b34f08ffSHong Zhang     if (ninfo > 40) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"number of INFO %d must <= 40\n",ninfo);
970b34f08ffSHong Zhang     ierr = PetscMalloc1(ninfo,&mumps->info);CHKERRQ(ierr);
971b34f08ffSHong Zhang     mumps->ninfo = ninfo;
972b34f08ffSHong Zhang     for (i=0; i<ninfo; i++) {
973b34f08ffSHong Zhang       if (info[i] < 0 || info[i]>40) {
974b34f08ffSHong Zhang         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index of INFO %d must between 1 and 40\n",ninfo);
975b34f08ffSHong Zhang       } else {
976b34f08ffSHong Zhang         mumps->info[i] = info[i];
977b34f08ffSHong Zhang       }
978b34f08ffSHong Zhang     }
979b34f08ffSHong Zhang   }
980b34f08ffSHong Zhang 
981dcd589f8SShri Abhyankar   PetscOptionsEnd();
982dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
983dcd589f8SShri Abhyankar }
984dcd589f8SShri Abhyankar 
985dcd589f8SShri Abhyankar #undef __FUNCT__
986dcd589f8SShri Abhyankar #define __FUNCT__ "PetscInitializeMUMPS"
987f697e70eSHong Zhang PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps)
988dcd589f8SShri Abhyankar {
989dcd589f8SShri Abhyankar   PetscErrorCode ierr;
990dcd589f8SShri Abhyankar 
991dcd589f8SShri Abhyankar   PetscFunctionBegin;
992ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)A), &mumps->myid);
993ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&mumps->size);CHKERRQ(ierr);
994ce94432eSBarry Smith   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(mumps->comm_mumps));CHKERRQ(ierr);
9952205254eSKarl Rupp 
996f697e70eSHong Zhang   mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps);
997f697e70eSHong Zhang 
998f697e70eSHong Zhang   mumps->id.job = JOB_INIT;
999f697e70eSHong Zhang   mumps->id.par = 1;  /* host participates factorizaton and solve */
1000f697e70eSHong Zhang   mumps->id.sym = mumps->sym;
10012907cef9SHong Zhang   PetscMUMPS_c(&mumps->id);
1002f697e70eSHong Zhang 
1003f697e70eSHong Zhang   mumps->CleanUpMUMPS = PETSC_FALSE;
10040298fd71SBarry Smith   mumps->scat_rhs     = NULL;
10050298fd71SBarry Smith   mumps->scat_sol     = NULL;
10069a2535b5SHong Zhang 
100770544d5fSHong Zhang   /* set PETSc-MUMPS default options - override MUMPS default */
10089a2535b5SHong Zhang   mumps->id.ICNTL(3) = 0;
10099a2535b5SHong Zhang   mumps->id.ICNTL(4) = 0;
10109a2535b5SHong Zhang   if (mumps->size == 1) {
10119a2535b5SHong Zhang     mumps->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
10129a2535b5SHong Zhang   } else {
10139a2535b5SHong Zhang     mumps->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
10144e34a73bSHong Zhang     mumps->id.ICNTL(20) = 0;   /* rhs is in dense format */
101570544d5fSHong Zhang     mumps->id.ICNTL(21) = 1;   /* distributed solution */
10169a2535b5SHong Zhang   }
1017dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
1018dcd589f8SShri Abhyankar }
1019dcd589f8SShri Abhyankar 
1020a5e57a09SHong 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 */
1021397b6df1SKris Buschelman #undef __FUNCT__
1022f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
10230481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1024b24902e0SBarry Smith {
1025a5e57a09SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
1026dcd589f8SShri Abhyankar   PetscErrorCode ierr;
102767877ebaSShri Abhyankar   Vec            b;
102867877ebaSShri Abhyankar   IS             is_iden;
102967877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1030397b6df1SKris Buschelman 
1031397b6df1SKris Buschelman   PetscFunctionBegin;
1032a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1033dcd589f8SShri Abhyankar 
10349a2535b5SHong Zhang   /* Set MUMPS options from the options database */
10359a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1036dcd589f8SShri Abhyankar 
1037a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1038dcd589f8SShri Abhyankar 
103967877ebaSShri Abhyankar   /* analysis phase */
104067877ebaSShri Abhyankar   /*----------------*/
1041a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1042a5e57a09SHong Zhang   mumps->id.n   = M;
1043a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
104467877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1045a5e57a09SHong Zhang     if (!mumps->myid) {
1046a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1047a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
1048*940cd9d6SSatish Balay         mumps->id.a = (MumpsScalar*)mumps->val;
104967877ebaSShri Abhyankar       }
1050a5e57a09SHong Zhang       if (mumps->id.ICNTL(7) == 1) { /* use user-provide matrix ordering - assuming r = c ordering */
10515248a706SHong Zhang         /*
10525248a706SHong Zhang         PetscBool      flag;
10535248a706SHong Zhang         ierr = ISEqual(r,c,&flag);CHKERRQ(ierr);
10545248a706SHong Zhang         if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"row_perm != col_perm");
10555248a706SHong Zhang         ierr = ISView(r,PETSC_VIEWER_STDOUT_SELF);
10565248a706SHong Zhang          */
1057a5e57a09SHong Zhang         if (!mumps->myid) {
1058e0b74bf9SHong Zhang           const PetscInt *idx;
1059e0b74bf9SHong Zhang           PetscInt       i,*perm_in;
10602205254eSKarl Rupp 
1061785e854fSJed Brown           ierr = PetscMalloc1(M,&perm_in);CHKERRQ(ierr);
1062e0b74bf9SHong Zhang           ierr = ISGetIndices(r,&idx);CHKERRQ(ierr);
10632205254eSKarl Rupp 
1064a5e57a09SHong Zhang           mumps->id.perm_in = perm_in;
1065e0b74bf9SHong Zhang           for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */
1066e0b74bf9SHong Zhang           ierr = ISRestoreIndices(r,&idx);CHKERRQ(ierr);
1067e0b74bf9SHong Zhang         }
1068e0b74bf9SHong Zhang       }
106967877ebaSShri Abhyankar     }
107067877ebaSShri Abhyankar     break;
107167877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1072a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1073a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1074a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
1075*940cd9d6SSatish Balay       mumps->id.a_loc = (MumpsScalar*)mumps->val;
107667877ebaSShri Abhyankar     }
107767877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1078a5e57a09SHong Zhang     if (!mumps->myid) {
10792cd7d884SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->rmap->N,&mumps->b_seq);CHKERRQ(ierr);
10802cd7d884SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->N,0,1,&is_iden);CHKERRQ(ierr);
108167877ebaSShri Abhyankar     } else {
1082a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
108367877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
108467877ebaSShri Abhyankar     }
10852a7a6963SBarry Smith     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1086a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
10876bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
10886bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
108967877ebaSShri Abhyankar     break;
109067877ebaSShri Abhyankar   }
1091a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
1092a5e57a09SHong 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));
109367877ebaSShri Abhyankar 
1094719d5645SBarry Smith   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1095dcd589f8SShri Abhyankar   F->ops->solve           = MatSolve_MUMPS;
109651d5961aSHong Zhang   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
10974e34a73bSHong Zhang   F->ops->matsolve        = MatMatSolve_MUMPS;
1098b24902e0SBarry Smith   PetscFunctionReturn(0);
1099b24902e0SBarry Smith }
1100b24902e0SBarry Smith 
1101450b117fSShri Abhyankar /* Note the Petsc r and c permutations are ignored */
1102450b117fSShri Abhyankar #undef __FUNCT__
1103450b117fSShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_BAIJMUMPS"
1104450b117fSShri Abhyankar PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1105450b117fSShri Abhyankar {
1106a5e57a09SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
1107dcd589f8SShri Abhyankar   PetscErrorCode ierr;
110867877ebaSShri Abhyankar   Vec            b;
110967877ebaSShri Abhyankar   IS             is_iden;
111067877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1111450b117fSShri Abhyankar 
1112450b117fSShri Abhyankar   PetscFunctionBegin;
1113a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1114dcd589f8SShri Abhyankar 
11159a2535b5SHong Zhang   /* Set MUMPS options from the options database */
11169a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1117dcd589f8SShri Abhyankar 
1118a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
111967877ebaSShri Abhyankar 
112067877ebaSShri Abhyankar   /* analysis phase */
112167877ebaSShri Abhyankar   /*----------------*/
1122a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1123a5e57a09SHong Zhang   mumps->id.n   = M;
1124a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
112567877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1126a5e57a09SHong Zhang     if (!mumps->myid) {
1127a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1128a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
1129*940cd9d6SSatish Balay         mumps->id.a = (MumpsScalar*)mumps->val;
113067877ebaSShri Abhyankar       }
113167877ebaSShri Abhyankar     }
113267877ebaSShri Abhyankar     break;
113367877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1134a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1135a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1136a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
1137*940cd9d6SSatish Balay       mumps->id.a_loc = (MumpsScalar*)mumps->val;
113867877ebaSShri Abhyankar     }
113967877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1140a5e57a09SHong Zhang     if (!mumps->myid) {
1141a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
114267877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
114367877ebaSShri Abhyankar     } else {
1144a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
114567877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
114667877ebaSShri Abhyankar     }
11472a7a6963SBarry Smith     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1148a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
11496bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
11506bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
115167877ebaSShri Abhyankar     break;
115267877ebaSShri Abhyankar   }
1153a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
1154a5e57a09SHong 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));
115567877ebaSShri Abhyankar 
1156450b117fSShri Abhyankar   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1157dcd589f8SShri Abhyankar   F->ops->solve           = MatSolve_MUMPS;
115851d5961aSHong Zhang   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1159450b117fSShri Abhyankar   PetscFunctionReturn(0);
1160450b117fSShri Abhyankar }
1161b24902e0SBarry Smith 
1162141f4205SHong Zhang /* Note the Petsc r permutation and factor info are ignored */
1163397b6df1SKris Buschelman #undef __FUNCT__
116467877ebaSShri Abhyankar #define __FUNCT__ "MatCholeskyFactorSymbolic_MUMPS"
116567877ebaSShri Abhyankar PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
1166b24902e0SBarry Smith {
1167a5e57a09SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
1168dcd589f8SShri Abhyankar   PetscErrorCode ierr;
116967877ebaSShri Abhyankar   Vec            b;
117067877ebaSShri Abhyankar   IS             is_iden;
117167877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1172397b6df1SKris Buschelman 
1173397b6df1SKris Buschelman   PetscFunctionBegin;
1174a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1175dcd589f8SShri Abhyankar 
11769a2535b5SHong Zhang   /* Set MUMPS options from the options database */
11779a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1178dcd589f8SShri Abhyankar 
1179a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1180dcd589f8SShri Abhyankar 
118167877ebaSShri Abhyankar   /* analysis phase */
118267877ebaSShri Abhyankar   /*----------------*/
1183a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1184a5e57a09SHong Zhang   mumps->id.n   = M;
1185a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
118667877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1187a5e57a09SHong Zhang     if (!mumps->myid) {
1188a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1189a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
1190*940cd9d6SSatish Balay         mumps->id.a = (MumpsScalar*)mumps->val;
119167877ebaSShri Abhyankar       }
119267877ebaSShri Abhyankar     }
119367877ebaSShri Abhyankar     break;
119467877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1195a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1196a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1197a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
1198*940cd9d6SSatish Balay       mumps->id.a_loc = (MumpsScalar*)mumps->val;
119967877ebaSShri Abhyankar     }
120067877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1201a5e57a09SHong Zhang     if (!mumps->myid) {
1202a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
120367877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
120467877ebaSShri Abhyankar     } else {
1205a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
120667877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
120767877ebaSShri Abhyankar     }
12082a7a6963SBarry Smith     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1209a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
12106bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
12116bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
121267877ebaSShri Abhyankar     break;
121367877ebaSShri Abhyankar   }
1214a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
1215a5e57a09SHong 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));
121667877ebaSShri Abhyankar 
12172792810eSHong Zhang   F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
1218dcd589f8SShri Abhyankar   F->ops->solve                 = MatSolve_MUMPS;
121951d5961aSHong Zhang   F->ops->solvetranspose        = MatSolve_MUMPS;
12204e34a73bSHong Zhang   F->ops->matsolve              = MatMatSolve_MUMPS;
12214e34a73bSHong Zhang #if defined(PETSC_USE_COMPLEX)
12220298fd71SBarry Smith   F->ops->getinertia = NULL;
12234e34a73bSHong Zhang #else
12244e34a73bSHong Zhang   F->ops->getinertia = MatGetInertia_SBAIJMUMPS;
1225db4efbfdSBarry Smith #endif
1226b24902e0SBarry Smith   PetscFunctionReturn(0);
1227b24902e0SBarry Smith }
1228b24902e0SBarry Smith 
1229397b6df1SKris Buschelman #undef __FUNCT__
123064e6c443SBarry Smith #define __FUNCT__ "MatView_MUMPS"
123164e6c443SBarry Smith PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
123274ed9c26SBarry Smith {
1233f6c57405SHong Zhang   PetscErrorCode    ierr;
123464e6c443SBarry Smith   PetscBool         iascii;
123564e6c443SBarry Smith   PetscViewerFormat format;
1236a5e57a09SHong Zhang   Mat_MUMPS         *mumps=(Mat_MUMPS*)A->spptr;
1237f6c57405SHong Zhang 
1238f6c57405SHong Zhang   PetscFunctionBegin;
123964e6c443SBarry Smith   /* check if matrix is mumps type */
124064e6c443SBarry Smith   if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0);
124164e6c443SBarry Smith 
1242251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
124364e6c443SBarry Smith   if (iascii) {
124464e6c443SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
124564e6c443SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO) {
124664e6c443SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
1247a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                   %d \n",mumps->id.sym);CHKERRQ(ierr);
1248a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):            %d \n",mumps->id.par);CHKERRQ(ierr);
1249a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):         %d \n",mumps->id.ICNTL(1));CHKERRQ(ierr);
1250a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg): %d \n",mumps->id.ICNTL(2));CHKERRQ(ierr);
1251a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):   %d \n",mumps->id.ICNTL(3));CHKERRQ(ierr);
1252a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):        %d \n",mumps->id.ICNTL(4));CHKERRQ(ierr);
1253a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):         %d \n",mumps->id.ICNTL(5));CHKERRQ(ierr);
1254a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):        %d \n",mumps->id.ICNTL(6));CHKERRQ(ierr);
1255a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (sequentia matrix ordering):%d \n",mumps->id.ICNTL(7));CHKERRQ(ierr);
1256a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scalling strategy):        %d \n",mumps->id.ICNTL(8));CHKERRQ(ierr);
1257a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements):  %d \n",mumps->id.ICNTL(10));CHKERRQ(ierr);
1258a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):          %d \n",mumps->id.ICNTL(11));CHKERRQ(ierr);
1259a5e57a09SHong Zhang       if (mumps->id.ICNTL(11)>0) {
1260a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(4) (inf norm of input mat):        %g\n",mumps->id.RINFOG(4));CHKERRQ(ierr);
1261a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(5) (inf norm of solution):         %g\n",mumps->id.RINFOG(5));CHKERRQ(ierr);
1262a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(6) (inf norm of residual):         %g\n",mumps->id.RINFOG(6));CHKERRQ(ierr);
1263a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8));CHKERRQ(ierr);
1264a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(9) (error estimate):               %g \n",mumps->id.RINFOG(9));CHKERRQ(ierr);
1265a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11));CHKERRQ(ierr);
1266f6c57405SHong Zhang       }
1267a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",mumps->id.ICNTL(12));CHKERRQ(ierr);
1268a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",mumps->id.ICNTL(13));CHKERRQ(ierr);
1269a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",mumps->id.ICNTL(14));CHKERRQ(ierr);
1270f6c57405SHong Zhang       /* ICNTL(15-17) not used */
1271a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",mumps->id.ICNTL(18));CHKERRQ(ierr);
1272a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Shur complement info):                       %d \n",mumps->id.ICNTL(19));CHKERRQ(ierr);
1273a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",mumps->id.ICNTL(20));CHKERRQ(ierr);
1274a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (somumpstion struct):                            %d \n",mumps->id.ICNTL(21));CHKERRQ(ierr);
1275a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",mumps->id.ICNTL(22));CHKERRQ(ierr);
1276a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",mumps->id.ICNTL(23));CHKERRQ(ierr);
1277c0165424SHong Zhang 
1278a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",mumps->id.ICNTL(24));CHKERRQ(ierr);
1279a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",mumps->id.ICNTL(25));CHKERRQ(ierr);
1280a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",mumps->id.ICNTL(26));CHKERRQ(ierr);
1281a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",mumps->id.ICNTL(27));CHKERRQ(ierr);
1282a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(28) (use parallel or sequential ordering):        %d \n",mumps->id.ICNTL(28));CHKERRQ(ierr);
1283a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(29) (parallel ordering):                          %d \n",mumps->id.ICNTL(29));CHKERRQ(ierr);
128442179a6aSHong Zhang 
1285a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(30) (user-specified set of entries in inv(A)):    %d \n",mumps->id.ICNTL(30));CHKERRQ(ierr);
1286a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(31) (factors is discarded in the solve phase):    %d \n",mumps->id.ICNTL(31));CHKERRQ(ierr);
1287a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(33) (compute determinant):                        %d \n",mumps->id.ICNTL(33));CHKERRQ(ierr);
1288f6c57405SHong Zhang 
1289a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",mumps->id.CNTL(1));CHKERRQ(ierr);
1290a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2));CHKERRQ(ierr);
1291a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absomumpste pivoting threshold):      %g \n",mumps->id.CNTL(3));CHKERRQ(ierr);
1292a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(4) (vamumpse of static pivoting):         %g \n",mumps->id.CNTL(4));CHKERRQ(ierr);
1293a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",mumps->id.CNTL(5));CHKERRQ(ierr);
1294f6c57405SHong Zhang 
1295f6c57405SHong Zhang       /* infomation local to each processor */
129634ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);
12977b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
1298a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %g \n",mumps->myid,mumps->id.RINFO(1));CHKERRQ(ierr);
129934ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
130034ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);
1301a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(2));CHKERRQ(ierr);
130234ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
130334ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);
1304a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(3));CHKERRQ(ierr);
130534ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1306f6c57405SHong Zhang 
130734ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr);
1308a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] %d \n",mumps->myid,mumps->id.INFO(15));CHKERRQ(ierr);
130934ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1310f6c57405SHong Zhang 
131134ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr);
1312a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(16));CHKERRQ(ierr);
131334ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1314f6c57405SHong Zhang 
131534ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);
1316a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(23));CHKERRQ(ierr);
131734ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1318b34f08ffSHong Zhang 
1319b34f08ffSHong Zhang       if (mumps->ninfo && mumps->ninfo <= 40){
1320b34f08ffSHong Zhang         PetscInt i;
1321b34f08ffSHong Zhang         for (i=0; i<mumps->ninfo; i++){
1322b34f08ffSHong Zhang           ierr = PetscViewerASCIIPrintf(viewer, "  INFO(%d): \n",mumps->info[i]);CHKERRQ(ierr);
1323b34f08ffSHong Zhang           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(mumps->info[i]));CHKERRQ(ierr);
1324b34f08ffSHong Zhang           ierr = PetscViewerFlush(viewer);
1325b34f08ffSHong Zhang         }
1326b34f08ffSHong Zhang       }
1327b34f08ffSHong Zhang 
1328b34f08ffSHong Zhang 
13297b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
1330f6c57405SHong Zhang 
1331a5e57a09SHong Zhang       if (!mumps->myid) { /* information from the host */
1332a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",mumps->id.RINFOG(1));CHKERRQ(ierr);
1333a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",mumps->id.RINFOG(2));CHKERRQ(ierr);
1334a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",mumps->id.RINFOG(3));CHKERRQ(ierr);
1335a5e57a09SHong 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);
1336f6c57405SHong Zhang 
1337a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(3));CHKERRQ(ierr);
1338a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(4));CHKERRQ(ierr);
1339a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",mumps->id.INFOG(5));CHKERRQ(ierr);
1340a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",mumps->id.INFOG(6));CHKERRQ(ierr);
1341a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively use after analysis): %d \n",mumps->id.INFOG(7));CHKERRQ(ierr);
1342a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",mumps->id.INFOG(8));CHKERRQ(ierr);
1343a5e57a09SHong 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);
1344a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",mumps->id.INFOG(10));CHKERRQ(ierr);
1345a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",mumps->id.INFOG(11));CHKERRQ(ierr);
1346a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",mumps->id.INFOG(12));CHKERRQ(ierr);
1347a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",mumps->id.INFOG(13));CHKERRQ(ierr);
1348a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",mumps->id.INFOG(14));CHKERRQ(ierr);
1349a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",mumps->id.INFOG(15));CHKERRQ(ierr);
1350a5e57a09SHong 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);
1351a5e57a09SHong 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);
1352a5e57a09SHong 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);
1353a5e57a09SHong 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);
1354a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",mumps->id.INFOG(20));CHKERRQ(ierr);
1355a5e57a09SHong 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);
1356a5e57a09SHong 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);
1357a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",mumps->id.INFOG(23));CHKERRQ(ierr);
1358a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",mumps->id.INFOG(24));CHKERRQ(ierr);
1359a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",mumps->id.INFOG(25));CHKERRQ(ierr);
136040d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(28) (after factorization: number of null pivots encountered): %d\n",mumps->id.INFOG(28));CHKERRQ(ierr);
136140d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(29) (after factorization: effective number of entries in the factors (sum over all processors)): %d\n",mumps->id.INFOG(29));CHKERRQ(ierr);
136240d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(30, 31) (after solution: size in Mbytes of memory used during solution phase): %d, %d\n",mumps->id.INFOG(30),mumps->id.INFOG(31));CHKERRQ(ierr);
136340d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(32) (after analysis: type of analysis done): %d\n",mumps->id.INFOG(32));CHKERRQ(ierr);
136440d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(33) (value used for ICNTL(8)): %d\n",mumps->id.INFOG(33));CHKERRQ(ierr);
136540d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(34) (exponent of the determinant if determinant is requested): %d\n",mumps->id.INFOG(34));CHKERRQ(ierr);
1366f6c57405SHong Zhang       }
1367f6c57405SHong Zhang     }
1368cb828f0fSHong Zhang   }
1369f6c57405SHong Zhang   PetscFunctionReturn(0);
1370f6c57405SHong Zhang }
1371f6c57405SHong Zhang 
137235bd34faSBarry Smith #undef __FUNCT__
137335bd34faSBarry Smith #define __FUNCT__ "MatGetInfo_MUMPS"
137435bd34faSBarry Smith PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
137535bd34faSBarry Smith {
1376cb828f0fSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)A->spptr;
137735bd34faSBarry Smith 
137835bd34faSBarry Smith   PetscFunctionBegin;
137935bd34faSBarry Smith   info->block_size        = 1.0;
1380cb828f0fSHong Zhang   info->nz_allocated      = mumps->id.INFOG(20);
1381cb828f0fSHong Zhang   info->nz_used           = mumps->id.INFOG(20);
138235bd34faSBarry Smith   info->nz_unneeded       = 0.0;
138335bd34faSBarry Smith   info->assemblies        = 0.0;
138435bd34faSBarry Smith   info->mallocs           = 0.0;
138535bd34faSBarry Smith   info->memory            = 0.0;
138635bd34faSBarry Smith   info->fill_ratio_given  = 0;
138735bd34faSBarry Smith   info->fill_ratio_needed = 0;
138835bd34faSBarry Smith   info->factor_mallocs    = 0;
138935bd34faSBarry Smith   PetscFunctionReturn(0);
139035bd34faSBarry Smith }
139135bd34faSBarry Smith 
13925ccb76cbSHong Zhang /* -------------------------------------------------------------------------------------------*/
13935ccb76cbSHong Zhang #undef __FUNCT__
13945ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl_MUMPS"
13955ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival)
13965ccb76cbSHong Zhang {
1397a5e57a09SHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
13985ccb76cbSHong Zhang 
13995ccb76cbSHong Zhang   PetscFunctionBegin;
1400a5e57a09SHong Zhang   mumps->id.ICNTL(icntl) = ival;
14015ccb76cbSHong Zhang   PetscFunctionReturn(0);
14025ccb76cbSHong Zhang }
14035ccb76cbSHong Zhang 
14045ccb76cbSHong Zhang #undef __FUNCT__
1405bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetIcntl_MUMPS"
1406bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt *ival)
1407bc6112feSHong Zhang {
1408bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1409bc6112feSHong Zhang 
1410bc6112feSHong Zhang   PetscFunctionBegin;
1411bc6112feSHong Zhang   *ival = mumps->id.ICNTL(icntl);
1412bc6112feSHong Zhang   PetscFunctionReturn(0);
1413bc6112feSHong Zhang }
1414bc6112feSHong Zhang 
1415bc6112feSHong Zhang #undef __FUNCT__
14165ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl"
14175ccb76cbSHong Zhang /*@
14185ccb76cbSHong Zhang   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()
14195ccb76cbSHong Zhang 
14205ccb76cbSHong Zhang    Logically Collective on Mat
14215ccb76cbSHong Zhang 
14225ccb76cbSHong Zhang    Input Parameters:
14235ccb76cbSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
14245ccb76cbSHong Zhang .  icntl - index of MUMPS parameter array ICNTL()
14255ccb76cbSHong Zhang -  ival - value of MUMPS ICNTL(icntl)
14265ccb76cbSHong Zhang 
14275ccb76cbSHong Zhang   Options Database:
14285ccb76cbSHong Zhang .   -mat_mumps_icntl_<icntl> <ival>
14295ccb76cbSHong Zhang 
14305ccb76cbSHong Zhang    Level: beginner
14315ccb76cbSHong Zhang 
14325ccb76cbSHong Zhang    References: MUMPS Users' Guide
14335ccb76cbSHong Zhang 
14345ccb76cbSHong Zhang .seealso: MatGetFactor()
14355ccb76cbSHong Zhang @*/
14365ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival)
14375ccb76cbSHong Zhang {
14385ccb76cbSHong Zhang   PetscErrorCode ierr;
14395ccb76cbSHong Zhang 
14405ccb76cbSHong Zhang   PetscFunctionBegin;
14415ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
14425ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,ival,3);
14435ccb76cbSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr);
14445ccb76cbSHong Zhang   PetscFunctionReturn(0);
14455ccb76cbSHong Zhang }
14465ccb76cbSHong Zhang 
1447bc6112feSHong Zhang #undef __FUNCT__
1448bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetIcntl"
1449a21f80fcSHong Zhang /*@
1450a21f80fcSHong Zhang   MatMumpsGetIcntl - Get MUMPS parameter ICNTL()
1451a21f80fcSHong Zhang 
1452a21f80fcSHong Zhang    Logically Collective on Mat
1453a21f80fcSHong Zhang 
1454a21f80fcSHong Zhang    Input Parameters:
1455a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1456a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array ICNTL()
1457a21f80fcSHong Zhang 
1458a21f80fcSHong Zhang   Output Parameter:
1459a21f80fcSHong Zhang .  ival - value of MUMPS ICNTL(icntl)
1460a21f80fcSHong Zhang 
1461a21f80fcSHong Zhang    Level: beginner
1462a21f80fcSHong Zhang 
1463a21f80fcSHong Zhang    References: MUMPS Users' Guide
1464a21f80fcSHong Zhang 
1465a21f80fcSHong Zhang .seealso: MatGetFactor()
1466a21f80fcSHong Zhang @*/
1467bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl(Mat F,PetscInt icntl,PetscInt *ival)
1468bc6112feSHong Zhang {
1469bc6112feSHong Zhang   PetscErrorCode ierr;
1470bc6112feSHong Zhang 
1471bc6112feSHong Zhang   PetscFunctionBegin;
1472bc6112feSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
1473bc6112feSHong Zhang   PetscValidIntPointer(ival,3);
1474bc6112feSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetIcntl_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
1475bc6112feSHong Zhang   PetscFunctionReturn(0);
1476bc6112feSHong Zhang }
1477bc6112feSHong Zhang 
14788928b65cSHong Zhang /* -------------------------------------------------------------------------------------------*/
14798928b65cSHong Zhang #undef __FUNCT__
14808928b65cSHong Zhang #define __FUNCT__ "MatMumpsSetCntl_MUMPS"
14818928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal val)
14828928b65cSHong Zhang {
14838928b65cSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
14848928b65cSHong Zhang 
14858928b65cSHong Zhang   PetscFunctionBegin;
14868928b65cSHong Zhang   mumps->id.CNTL(icntl) = val;
14878928b65cSHong Zhang   PetscFunctionReturn(0);
14888928b65cSHong Zhang }
14898928b65cSHong Zhang 
14908928b65cSHong Zhang #undef __FUNCT__
1491bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetCntl_MUMPS"
1492bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal *val)
1493bc6112feSHong Zhang {
1494bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1495bc6112feSHong Zhang 
1496bc6112feSHong Zhang   PetscFunctionBegin;
1497bc6112feSHong Zhang   *val = mumps->id.CNTL(icntl);
1498bc6112feSHong Zhang   PetscFunctionReturn(0);
1499bc6112feSHong Zhang }
1500bc6112feSHong Zhang 
1501bc6112feSHong Zhang #undef __FUNCT__
15028928b65cSHong Zhang #define __FUNCT__ "MatMumpsSetCntl"
15038928b65cSHong Zhang /*@
15048928b65cSHong Zhang   MatMumpsSetCntl - Set MUMPS parameter CNTL()
15058928b65cSHong Zhang 
15068928b65cSHong Zhang    Logically Collective on Mat
15078928b65cSHong Zhang 
15088928b65cSHong Zhang    Input Parameters:
15098928b65cSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
15108928b65cSHong Zhang .  icntl - index of MUMPS parameter array CNTL()
15118928b65cSHong Zhang -  val - value of MUMPS CNTL(icntl)
15128928b65cSHong Zhang 
15138928b65cSHong Zhang   Options Database:
15148928b65cSHong Zhang .   -mat_mumps_cntl_<icntl> <val>
15158928b65cSHong Zhang 
15168928b65cSHong Zhang    Level: beginner
15178928b65cSHong Zhang 
15188928b65cSHong Zhang    References: MUMPS Users' Guide
15198928b65cSHong Zhang 
15208928b65cSHong Zhang .seealso: MatGetFactor()
15218928b65cSHong Zhang @*/
15228928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl(Mat F,PetscInt icntl,PetscReal val)
15238928b65cSHong Zhang {
15248928b65cSHong Zhang   PetscErrorCode ierr;
15258928b65cSHong Zhang 
15268928b65cSHong Zhang   PetscFunctionBegin;
15278928b65cSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
1528bc6112feSHong Zhang   PetscValidLogicalCollectiveReal(F,val,3);
15298928b65cSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsSetCntl_C",(Mat,PetscInt,PetscReal),(F,icntl,val));CHKERRQ(ierr);
15308928b65cSHong Zhang   PetscFunctionReturn(0);
15318928b65cSHong Zhang }
15328928b65cSHong Zhang 
1533bc6112feSHong Zhang #undef __FUNCT__
1534bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetCntl"
1535a21f80fcSHong Zhang /*@
1536a21f80fcSHong Zhang   MatMumpsGetCntl - Get MUMPS parameter CNTL()
1537a21f80fcSHong Zhang 
1538a21f80fcSHong Zhang    Logically Collective on Mat
1539a21f80fcSHong Zhang 
1540a21f80fcSHong Zhang    Input Parameters:
1541a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1542a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array CNTL()
1543a21f80fcSHong Zhang 
1544a21f80fcSHong Zhang   Output Parameter:
1545a21f80fcSHong Zhang .  val - value of MUMPS CNTL(icntl)
1546a21f80fcSHong Zhang 
1547a21f80fcSHong Zhang    Level: beginner
1548a21f80fcSHong Zhang 
1549a21f80fcSHong Zhang    References: MUMPS Users' Guide
1550a21f80fcSHong Zhang 
1551a21f80fcSHong Zhang .seealso: MatGetFactor()
1552a21f80fcSHong Zhang @*/
1553bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl(Mat F,PetscInt icntl,PetscReal *val)
1554bc6112feSHong Zhang {
1555bc6112feSHong Zhang   PetscErrorCode ierr;
1556bc6112feSHong Zhang 
1557bc6112feSHong Zhang   PetscFunctionBegin;
1558bc6112feSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
1559bc6112feSHong Zhang   PetscValidRealPointer(val,3);
1560bc6112feSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetCntl_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
1561bc6112feSHong Zhang   PetscFunctionReturn(0);
1562bc6112feSHong Zhang }
1563bc6112feSHong Zhang 
1564bc6112feSHong Zhang #undef __FUNCT__
1565ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfo_MUMPS"
1566ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F,PetscInt icntl,PetscInt *info)
1567bc6112feSHong Zhang {
1568bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1569bc6112feSHong Zhang 
1570bc6112feSHong Zhang   PetscFunctionBegin;
1571bc6112feSHong Zhang   *info = mumps->id.INFO(icntl);
1572bc6112feSHong Zhang   PetscFunctionReturn(0);
1573bc6112feSHong Zhang }
1574bc6112feSHong Zhang 
1575bc6112feSHong Zhang #undef __FUNCT__
1576ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfog_MUMPS"
1577ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F,PetscInt icntl,PetscInt *infog)
1578bc6112feSHong Zhang {
1579bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1580bc6112feSHong Zhang 
1581bc6112feSHong Zhang   PetscFunctionBegin;
1582bc6112feSHong Zhang   *infog = mumps->id.INFOG(icntl);
1583bc6112feSHong Zhang   PetscFunctionReturn(0);
1584bc6112feSHong Zhang }
1585bc6112feSHong Zhang 
1586bc6112feSHong Zhang #undef __FUNCT__
1587ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfo_MUMPS"
1588ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfo)
1589bc6112feSHong Zhang {
1590bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1591bc6112feSHong Zhang 
1592bc6112feSHong Zhang   PetscFunctionBegin;
1593bc6112feSHong Zhang   *rinfo = mumps->id.RINFO(icntl);
1594bc6112feSHong Zhang   PetscFunctionReturn(0);
1595bc6112feSHong Zhang }
1596bc6112feSHong Zhang 
1597bc6112feSHong Zhang #undef __FUNCT__
1598ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfog_MUMPS"
1599ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfog)
1600bc6112feSHong Zhang {
1601bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1602bc6112feSHong Zhang 
1603bc6112feSHong Zhang   PetscFunctionBegin;
1604bc6112feSHong Zhang   *rinfog = mumps->id.RINFOG(icntl);
1605bc6112feSHong Zhang   PetscFunctionReturn(0);
1606bc6112feSHong Zhang }
1607bc6112feSHong Zhang 
1608bc6112feSHong Zhang #undef __FUNCT__
1609ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfo"
1610a21f80fcSHong Zhang /*@
1611a21f80fcSHong Zhang   MatMumpsGetInfo - Get MUMPS parameter INFO()
1612a21f80fcSHong Zhang 
1613a21f80fcSHong Zhang    Logically Collective on Mat
1614a21f80fcSHong Zhang 
1615a21f80fcSHong Zhang    Input Parameters:
1616a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1617a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array INFO()
1618a21f80fcSHong Zhang 
1619a21f80fcSHong Zhang   Output Parameter:
1620a21f80fcSHong Zhang .  ival - value of MUMPS INFO(icntl)
1621a21f80fcSHong Zhang 
1622a21f80fcSHong Zhang    Level: beginner
1623a21f80fcSHong Zhang 
1624a21f80fcSHong Zhang    References: MUMPS Users' Guide
1625a21f80fcSHong Zhang 
1626a21f80fcSHong Zhang .seealso: MatGetFactor()
1627a21f80fcSHong Zhang @*/
1628ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo(Mat F,PetscInt icntl,PetscInt *ival)
1629bc6112feSHong Zhang {
1630bc6112feSHong Zhang   PetscErrorCode ierr;
1631bc6112feSHong Zhang 
1632bc6112feSHong Zhang   PetscFunctionBegin;
1633ca810319SHong Zhang   PetscValidIntPointer(ival,3);
1634ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetInfo_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
1635bc6112feSHong Zhang   PetscFunctionReturn(0);
1636bc6112feSHong Zhang }
1637bc6112feSHong Zhang 
1638bc6112feSHong Zhang #undef __FUNCT__
1639ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfog"
1640a21f80fcSHong Zhang /*@
1641a21f80fcSHong Zhang   MatMumpsGetInfog - Get MUMPS parameter INFOG()
1642a21f80fcSHong Zhang 
1643a21f80fcSHong Zhang    Logically Collective on Mat
1644a21f80fcSHong Zhang 
1645a21f80fcSHong Zhang    Input Parameters:
1646a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1647a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array INFOG()
1648a21f80fcSHong Zhang 
1649a21f80fcSHong Zhang   Output Parameter:
1650a21f80fcSHong Zhang .  ival - value of MUMPS INFOG(icntl)
1651a21f80fcSHong Zhang 
1652a21f80fcSHong Zhang    Level: beginner
1653a21f80fcSHong Zhang 
1654a21f80fcSHong Zhang    References: MUMPS Users' Guide
1655a21f80fcSHong Zhang 
1656a21f80fcSHong Zhang .seealso: MatGetFactor()
1657a21f80fcSHong Zhang @*/
1658ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog(Mat F,PetscInt icntl,PetscInt *ival)
1659bc6112feSHong Zhang {
1660bc6112feSHong Zhang   PetscErrorCode ierr;
1661bc6112feSHong Zhang 
1662bc6112feSHong Zhang   PetscFunctionBegin;
1663ca810319SHong Zhang   PetscValidIntPointer(ival,3);
1664ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetInfog_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
1665bc6112feSHong Zhang   PetscFunctionReturn(0);
1666bc6112feSHong Zhang }
1667bc6112feSHong Zhang 
1668bc6112feSHong Zhang #undef __FUNCT__
1669ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfo"
1670a21f80fcSHong Zhang /*@
1671a21f80fcSHong Zhang   MatMumpsGetRinfo - Get MUMPS parameter RINFO()
1672a21f80fcSHong Zhang 
1673a21f80fcSHong Zhang    Logically Collective on Mat
1674a21f80fcSHong Zhang 
1675a21f80fcSHong Zhang    Input Parameters:
1676a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1677a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array RINFO()
1678a21f80fcSHong Zhang 
1679a21f80fcSHong Zhang   Output Parameter:
1680a21f80fcSHong Zhang .  val - value of MUMPS RINFO(icntl)
1681a21f80fcSHong Zhang 
1682a21f80fcSHong Zhang    Level: beginner
1683a21f80fcSHong Zhang 
1684a21f80fcSHong Zhang    References: MUMPS Users' Guide
1685a21f80fcSHong Zhang 
1686a21f80fcSHong Zhang .seealso: MatGetFactor()
1687a21f80fcSHong Zhang @*/
1688ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo(Mat F,PetscInt icntl,PetscReal *val)
1689bc6112feSHong Zhang {
1690bc6112feSHong Zhang   PetscErrorCode ierr;
1691bc6112feSHong Zhang 
1692bc6112feSHong Zhang   PetscFunctionBegin;
1693bc6112feSHong Zhang   PetscValidRealPointer(val,3);
1694ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetRinfo_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
1695bc6112feSHong Zhang   PetscFunctionReturn(0);
1696bc6112feSHong Zhang }
1697bc6112feSHong Zhang 
1698bc6112feSHong Zhang #undef __FUNCT__
1699ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfog"
1700a21f80fcSHong Zhang /*@
1701a21f80fcSHong Zhang   MatMumpsGetRinfog - Get MUMPS parameter RINFOG()
1702a21f80fcSHong Zhang 
1703a21f80fcSHong Zhang    Logically Collective on Mat
1704a21f80fcSHong Zhang 
1705a21f80fcSHong Zhang    Input Parameters:
1706a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1707a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array RINFOG()
1708a21f80fcSHong Zhang 
1709a21f80fcSHong Zhang   Output Parameter:
1710a21f80fcSHong Zhang .  val - value of MUMPS RINFOG(icntl)
1711a21f80fcSHong Zhang 
1712a21f80fcSHong Zhang    Level: beginner
1713a21f80fcSHong Zhang 
1714a21f80fcSHong Zhang    References: MUMPS Users' Guide
1715a21f80fcSHong Zhang 
1716a21f80fcSHong Zhang .seealso: MatGetFactor()
1717a21f80fcSHong Zhang @*/
1718ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog(Mat F,PetscInt icntl,PetscReal *val)
1719bc6112feSHong Zhang {
1720bc6112feSHong Zhang   PetscErrorCode ierr;
1721bc6112feSHong Zhang 
1722bc6112feSHong Zhang   PetscFunctionBegin;
1723bc6112feSHong Zhang   PetscValidRealPointer(val,3);
1724ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetRinfog_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
1725bc6112feSHong Zhang   PetscFunctionReturn(0);
1726bc6112feSHong Zhang }
1727bc6112feSHong Zhang 
172824b6179bSKris Buschelman /*MC
17292692d6eeSBarry Smith   MATSOLVERMUMPS -  A matrix type providing direct solvers (LU and Cholesky) for
173024b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
173124b6179bSKris Buschelman 
173241c8de11SBarry Smith   Works with MATAIJ and MATSBAIJ matrices
173324b6179bSKris Buschelman 
173424b6179bSKris Buschelman   Options Database Keys:
17354e34a73bSHong Zhang +  -mat_mumps_icntl_1 <6>: ICNTL(1): output stream for error messages (None)
17364e34a73bSHong Zhang .  -mat_mumps_icntl_2 <0>: ICNTL(2): output stream for diagnostic printing, statistics, and warning (None)
17374e34a73bSHong Zhang .  -mat_mumps_icntl_3 <0>: ICNTL(3): output stream for global information, collected on the host (None)
17384e34a73bSHong Zhang .  -mat_mumps_icntl_4 <0>: ICNTL(4): level of printing (0 to 4) (None)
17394e34a73bSHong Zhang .  -mat_mumps_icntl_6 <7>: ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7) (None)
17404e34a73bSHong Zhang .  -mat_mumps_icntl_7 <7>: ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis (None)
17414e34a73bSHong Zhang .  -mat_mumps_icntl_8 <77>: ICNTL(8): scaling strategy (-2 to 8 or 77) (None)
17424e34a73bSHong Zhang .  -mat_mumps_icntl_10 <0>: ICNTL(10): max num of refinements (None)
17434e34a73bSHong Zhang .  -mat_mumps_icntl_11 <0>: ICNTL(11): statistics related to an error analysis (via -ksp_view) (None)
17444e34a73bSHong Zhang .  -mat_mumps_icntl_12 <1>: ICNTL(12): an ordering strategy for symmetric matrices (0 to 3) (None)
17454e34a73bSHong Zhang .  -mat_mumps_icntl_13 <0>: ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting (None)
17464e34a73bSHong Zhang .  -mat_mumps_icntl_14 <20>: ICNTL(14): percentage increase in the estimated working space (None)
17474e34a73bSHong Zhang .  -mat_mumps_icntl_19 <0>: ICNTL(19): computes the Schur complement (None)
17484e34a73bSHong Zhang .  -mat_mumps_icntl_22 <0>: ICNTL(22): in-core/out-of-core factorization and solve (0 or 1) (None)
17494e34a73bSHong Zhang .  -mat_mumps_icntl_23 <0>: ICNTL(23): max size of the working memory (MB) that can allocate per processor (None)
17504e34a73bSHong Zhang .  -mat_mumps_icntl_24 <0>: ICNTL(24): detection of null pivot rows (0 or 1) (None)
17514e34a73bSHong Zhang .  -mat_mumps_icntl_25 <0>: ICNTL(25): compute a solution of a deficient matrix and a null space basis (None)
17524e34a73bSHong Zhang .  -mat_mumps_icntl_26 <0>: ICNTL(26): drives the solution phase if a Schur complement matrix (None)
17534e34a73bSHong Zhang .  -mat_mumps_icntl_28 <1>: ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering (None)
17544e34a73bSHong Zhang .  -mat_mumps_icntl_29 <0>: ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis (None)
17554e34a73bSHong Zhang .  -mat_mumps_icntl_30 <0>: ICNTL(30): compute user-specified set of entries in inv(A) (None)
17564e34a73bSHong Zhang .  -mat_mumps_icntl_31 <0>: ICNTL(31): indicates which factors may be discarded during factorization (None)
17574e34a73bSHong Zhang .  -mat_mumps_icntl_33 <0>: ICNTL(33): compute determinant (None)
17584e34a73bSHong Zhang .  -mat_mumps_cntl_1 <0.01>: CNTL(1): relative pivoting threshold (None)
17594e34a73bSHong Zhang .  -mat_mumps_cntl_2 <1.49012e-08>: CNTL(2): stopping criterion of refinement (None)
17604e34a73bSHong Zhang .  -mat_mumps_cntl_3 <0>: CNTL(3): absolute pivoting threshold (None)
17614e34a73bSHong Zhang .  -mat_mumps_cntl_4 <-1>: CNTL(4): value for static pivoting (None)
17624e34a73bSHong Zhang -  -mat_mumps_cntl_5 <0>: CNTL(5): fixation for null pivots (None)
176324b6179bSKris Buschelman 
176424b6179bSKris Buschelman   Level: beginner
176524b6179bSKris Buschelman 
176641c8de11SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage
176741c8de11SBarry Smith 
176824b6179bSKris Buschelman M*/
176924b6179bSKris Buschelman 
177035bd34faSBarry Smith #undef __FUNCT__
177135bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_mumps"
1772f7a08781SBarry Smith static PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type)
177335bd34faSBarry Smith {
177435bd34faSBarry Smith   PetscFunctionBegin;
17752692d6eeSBarry Smith   *type = MATSOLVERMUMPS;
177635bd34faSBarry Smith   PetscFunctionReturn(0);
177735bd34faSBarry Smith }
177835bd34faSBarry Smith 
1779bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI AIJ matrices */
17802877fffaSHong Zhang #undef __FUNCT__
1781bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_aij_mumps"
17828cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
17832877fffaSHong Zhang {
17842877fffaSHong Zhang   Mat            B;
17852877fffaSHong Zhang   PetscErrorCode ierr;
17862877fffaSHong Zhang   Mat_MUMPS      *mumps;
1787ace3abfcSBarry Smith   PetscBool      isSeqAIJ;
17882877fffaSHong Zhang 
17892877fffaSHong Zhang   PetscFunctionBegin;
17902877fffaSHong Zhang   /* Create the factorization matrix */
1791251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr);
1792ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
17932877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
17942877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1795bccb9932SShri Abhyankar   if (isSeqAIJ) {
17960298fd71SBarry Smith     ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr);
1797bccb9932SShri Abhyankar   } else {
17980298fd71SBarry Smith     ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr);
1799bccb9932SShri Abhyankar   }
18002877fffaSHong Zhang 
1801b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
18022205254eSKarl Rupp 
18032877fffaSHong Zhang   B->ops->view        = MatView_MUMPS;
180435bd34faSBarry Smith   B->ops->getinfo     = MatGetInfo_MUMPS;
180520be8e61SHong Zhang   B->ops->getdiagonal = MatGetDiagonal_MUMPS;
18062205254eSKarl Rupp 
1807bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
1808bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
1809bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
1810bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
1811bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
1812bc6112feSHong Zhang 
1813ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
1814ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
1815ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
1816ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
1817450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
1818450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
1819d5f3da31SBarry Smith     B->factortype            = MAT_FACTOR_LU;
1820bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
1821bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
1822746480a1SHong Zhang     mumps->sym = 0;
1823dcd589f8SShri Abhyankar   } else {
182467877ebaSShri Abhyankar     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
1825450b117fSShri Abhyankar     B->factortype                  = MAT_FACTOR_CHOLESKY;
1826bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
1827bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
18286fdc2a6dSBarry Smith     if (A->spd_set && A->spd) mumps->sym = 1;
18296fdc2a6dSBarry Smith     else                      mumps->sym = 2;
1830450b117fSShri Abhyankar   }
18312877fffaSHong Zhang 
18322877fffaSHong Zhang   mumps->isAIJ    = PETSC_TRUE;
1833bf0cc555SLisandro Dalcin   mumps->Destroy  = B->ops->destroy;
18342877fffaSHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
18352877fffaSHong Zhang   B->spptr        = (void*)mumps;
18362205254eSKarl Rupp 
1837f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1838746480a1SHong Zhang 
18392877fffaSHong Zhang   *F = B;
18402877fffaSHong Zhang   PetscFunctionReturn(0);
18412877fffaSHong Zhang }
18422877fffaSHong Zhang 
1843bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI SBAIJ matrices */
18442877fffaSHong Zhang #undef __FUNCT__
1845bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_sbaij_mumps"
18468cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
18472877fffaSHong Zhang {
18482877fffaSHong Zhang   Mat            B;
18492877fffaSHong Zhang   PetscErrorCode ierr;
18502877fffaSHong Zhang   Mat_MUMPS      *mumps;
1851ace3abfcSBarry Smith   PetscBool      isSeqSBAIJ;
18522877fffaSHong Zhang 
18532877fffaSHong Zhang   PetscFunctionBegin;
1854ce94432eSBarry Smith   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
1855ce94432eSBarry Smith   if (A->rmap->bs > 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with block size > 1 with MUMPS Cholesky, use AIJ matrix instead");
1856251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
18572877fffaSHong Zhang   /* Create the factorization matrix */
1858ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
18592877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
18602877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1861b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
1862bccb9932SShri Abhyankar   if (isSeqSBAIJ) {
18630298fd71SBarry Smith     ierr = MatSeqSBAIJSetPreallocation(B,1,0,NULL);CHKERRQ(ierr);
18642205254eSKarl Rupp 
186516ebf90aSShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
1866dcd589f8SShri Abhyankar   } else {
18670298fd71SBarry Smith     ierr = MatMPISBAIJSetPreallocation(B,1,0,NULL,0,NULL);CHKERRQ(ierr);
18682205254eSKarl Rupp 
1869bccb9932SShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
1870bccb9932SShri Abhyankar   }
1871bccb9932SShri Abhyankar 
187267877ebaSShri Abhyankar   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
1873bccb9932SShri Abhyankar   B->ops->view                   = MatView_MUMPS;
187420be8e61SHong Zhang   B->ops->getdiagonal            = MatGetDiagonal_MUMPS;
18752205254eSKarl Rupp 
1876bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
1877b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
1878b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
1879b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
1880b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
1881bc6112feSHong Zhang 
1882ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
1883ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
1884ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
1885ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
18862205254eSKarl Rupp 
1887f4762488SHong Zhang   B->factortype = MAT_FACTOR_CHOLESKY;
18886fdc2a6dSBarry Smith   if (A->spd_set && A->spd) mumps->sym = 1;
18896fdc2a6dSBarry Smith   else                      mumps->sym = 2;
1890a214ac2aSShri Abhyankar 
1891bccb9932SShri Abhyankar   mumps->isAIJ    = PETSC_FALSE;
1892bf0cc555SLisandro Dalcin   mumps->Destroy  = B->ops->destroy;
1893f3c0ef26SHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
18942877fffaSHong Zhang   B->spptr        = (void*)mumps;
18952205254eSKarl Rupp 
1896f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1897746480a1SHong Zhang 
18982877fffaSHong Zhang   *F = B;
18992877fffaSHong Zhang   PetscFunctionReturn(0);
19002877fffaSHong Zhang }
190197969023SHong Zhang 
1902450b117fSShri Abhyankar #undef __FUNCT__
1903bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_baij_mumps"
19048cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
190567877ebaSShri Abhyankar {
190667877ebaSShri Abhyankar   Mat            B;
190767877ebaSShri Abhyankar   PetscErrorCode ierr;
190867877ebaSShri Abhyankar   Mat_MUMPS      *mumps;
1909ace3abfcSBarry Smith   PetscBool      isSeqBAIJ;
191067877ebaSShri Abhyankar 
191167877ebaSShri Abhyankar   PetscFunctionBegin;
191267877ebaSShri Abhyankar   /* Create the factorization matrix */
1913251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr);
1914ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
191567877ebaSShri Abhyankar   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
191667877ebaSShri Abhyankar   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1917bccb9932SShri Abhyankar   if (isSeqBAIJ) {
19180298fd71SBarry Smith     ierr = MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,NULL);CHKERRQ(ierr);
1919bccb9932SShri Abhyankar   } else {
19200298fd71SBarry Smith     ierr = MatMPIBAIJSetPreallocation(B,A->rmap->bs,0,NULL,0,NULL);CHKERRQ(ierr);
1921bccb9932SShri Abhyankar   }
1922450b117fSShri Abhyankar 
1923b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
1924450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
1925450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
1926450b117fSShri Abhyankar     B->factortype            = MAT_FACTOR_LU;
1927bccb9932SShri Abhyankar     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
1928bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
1929746480a1SHong Zhang     mumps->sym = 0;
1930f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n");
1931bccb9932SShri Abhyankar 
1932450b117fSShri Abhyankar   B->ops->view        = MatView_MUMPS;
193320be8e61SHong Zhang   B->ops->getdiagonal = MatGetDiagonal_MUMPS;
19342205254eSKarl Rupp 
1935bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
1936bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
1937bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
1938bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
1939bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
1940bc6112feSHong Zhang 
1941ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
1942ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
1943ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
1944ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
1945450b117fSShri Abhyankar 
1946450b117fSShri Abhyankar   mumps->isAIJ    = PETSC_TRUE;
1947bf0cc555SLisandro Dalcin   mumps->Destroy  = B->ops->destroy;
1948450b117fSShri Abhyankar   B->ops->destroy = MatDestroy_MUMPS;
1949450b117fSShri Abhyankar   B->spptr        = (void*)mumps;
19502205254eSKarl Rupp 
1951f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1952746480a1SHong Zhang 
1953450b117fSShri Abhyankar   *F = B;
1954450b117fSShri Abhyankar   PetscFunctionReturn(0);
1955450b117fSShri Abhyankar }
195642c9c57cSBarry Smith 
195742c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mumps(Mat,MatFactorType,Mat*);
195842c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_baij_mumps(Mat,MatFactorType,Mat*);
195942c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat,MatFactorType,Mat*);
196042c9c57cSBarry Smith 
196142c9c57cSBarry Smith #undef __FUNCT__
196242c9c57cSBarry Smith #define __FUNCT__ "MatSolverPackageRegister_MUMPS"
196329b38603SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_MUMPS(void)
196442c9c57cSBarry Smith {
196542c9c57cSBarry Smith   PetscErrorCode ierr;
196642c9c57cSBarry Smith 
196742c9c57cSBarry Smith   PetscFunctionBegin;
196842c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ,         MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr);
196942c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ,         MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr);
197042c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ,         MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr);
197142c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ,         MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr);
197242c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPISBAIJ,         MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr);
197342c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ,         MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr);
197442c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ,         MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr);
197542c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ,         MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr);
197642c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ,         MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr);
197742c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQSBAIJ,         MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr);
197842c9c57cSBarry Smith   PetscFunctionReturn(0);
197942c9c57cSBarry Smith }
198042c9c57cSBarry Smith 
1981