xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision 2cd7d8848ee56197cb53f945aadfa2185658b8d6)
11c2a3de1SBarry Smith 
2397b6df1SKris Buschelman /*
3c2b5dc30SHong Zhang     Provides an interface to the MUMPS sparse solver
4397b6df1SKris Buschelman */
551d5961aSHong Zhang 
6c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> /*I  "petscmat.h"  I*/
7c6db04a5SJed Brown #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
8397b6df1SKris Buschelman 
9397b6df1SKris Buschelman EXTERN_C_BEGIN
10397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
112907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
122907cef9SHong Zhang #include <cmumps_c.h>
132907cef9SHong Zhang #else
14c6db04a5SJed Brown #include <zmumps_c.h>
152907cef9SHong Zhang #endif
162907cef9SHong Zhang #else
172907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
182907cef9SHong Zhang #include <smumps_c.h>
19397b6df1SKris Buschelman #else
20c6db04a5SJed Brown #include <dmumps_c.h>
21397b6df1SKris Buschelman #endif
222907cef9SHong Zhang #endif
23397b6df1SKris Buschelman EXTERN_C_END
24397b6df1SKris Buschelman #define JOB_INIT -1
253d472b54SHong Zhang #define JOB_FACTSYMBOLIC 1
263d472b54SHong Zhang #define JOB_FACTNUMERIC 2
273d472b54SHong Zhang #define JOB_SOLVE 3
28397b6df1SKris Buschelman #define JOB_END -2
293d472b54SHong Zhang 
302907cef9SHong Zhang /* calls to MUMPS */
312907cef9SHong Zhang #if defined(PETSC_USE_COMPLEX)
322907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
332907cef9SHong Zhang #define PetscMUMPS_c cmumps_c
342907cef9SHong Zhang #else
352907cef9SHong Zhang #define PetscMUMPS_c zmumps_c
362907cef9SHong Zhang #endif
372907cef9SHong Zhang #else
382907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
392907cef9SHong Zhang #define PetscMUMPS_c smumps_c
402907cef9SHong Zhang #else
412907cef9SHong Zhang #define PetscMUMPS_c dmumps_c
422907cef9SHong Zhang #endif
432907cef9SHong Zhang #endif
442907cef9SHong Zhang 
453d472b54SHong Zhang 
46397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */
47397b6df1SKris Buschelman #define ICNTL(I) icntl[(I)-1]
48397b6df1SKris Buschelman #define CNTL(I) cntl[(I)-1]
49397b6df1SKris Buschelman #define INFOG(I) infog[(I)-1]
50a7aca84bSHong Zhang #define INFO(I) info[(I)-1]
51397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1]
52adc1d99fSHong Zhang #define RINFO(I) rinfo[(I)-1]
53397b6df1SKris Buschelman 
54397b6df1SKris Buschelman typedef struct {
55397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
562907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
572907cef9SHong Zhang   CMUMPS_STRUC_C id;
582907cef9SHong Zhang #else
59397b6df1SKris Buschelman   ZMUMPS_STRUC_C id;
602907cef9SHong Zhang #endif
612907cef9SHong Zhang #else
622907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
632907cef9SHong Zhang   SMUMPS_STRUC_C id;
64397b6df1SKris Buschelman #else
65397b6df1SKris Buschelman   DMUMPS_STRUC_C id;
66397b6df1SKris Buschelman #endif
672907cef9SHong Zhang #endif
682907cef9SHong Zhang 
69397b6df1SKris Buschelman   MatStructure matstruc;
70c1490034SHong Zhang   PetscMPIInt  myid,size;
71a5e57a09SHong Zhang   PetscInt     *irn,*jcn,nz,sym;
72397b6df1SKris Buschelman   PetscScalar  *val;
73397b6df1SKris Buschelman   MPI_Comm     comm_mumps;
74329ec9b3SHong Zhang   VecScatter   scat_rhs, scat_sol;
7564e6c443SBarry Smith   PetscBool    isAIJ,CleanUpMUMPS;
76*2cd7d884SHong Zhang   Vec          b_seq,x_seq;          /* used by MatSolve() */
77*2cd7d884SHong Zhang   Vec          bb_mpi,bb_seq,xx_mpi; /* used by MatMatSolve() */
78a5e57a09SHong Zhang   PetscInt     ICNTL9_pre;   /* check if ICNTL(9) is changed from previous MatSolve */
792205254eSKarl Rupp 
80bf0cc555SLisandro Dalcin   PetscErrorCode (*Destroy)(Mat);
81bccb9932SShri Abhyankar   PetscErrorCode (*ConvertToTriples)(Mat, int, MatReuse, int*, int**, int**, PetscScalar**);
82f0c56d0fSKris Buschelman } Mat_MUMPS;
83f0c56d0fSKris Buschelman 
8409573ac7SBarry Smith extern PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);
85b24902e0SBarry Smith 
86397b6df1SKris Buschelman /*
87d341cd04SHong Zhang   MatConvertToTriples_A_B - convert Petsc matrix to triples: row[nz], col[nz], val[nz]
88d341cd04SHong Zhang 
89397b6df1SKris Buschelman   input:
9067877ebaSShri Abhyankar     A       - matrix in aij,baij or sbaij (bs=1) format
91397b6df1SKris Buschelman     shift   - 0: C style output triple; 1: Fortran style output triple.
92bccb9932SShri Abhyankar     reuse   - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple
93bccb9932SShri Abhyankar               MAT_REUSE_MATRIX:   only the values in v array are updated
94397b6df1SKris Buschelman   output:
95397b6df1SKris Buschelman     nnz     - dim of r, c, and v (number of local nonzero entries of A)
96397b6df1SKris Buschelman     r, c, v - row and col index, matrix values (matrix triples)
97eb9baa12SBarry Smith 
98eb9baa12SBarry Smith   The returned values r, c, and sometimes v are obtained in a single PetscMalloc(). Then in MatDestroy_MUMPS() it is
99eb9baa12SBarry Smith   freed with PetscFree((mumps->irn);  This is not ideal code, the fact that v is ONLY sometimes part of mumps->irn means
100eb9baa12SBarry Smith   that the PetscMalloc() cannot easily be replaced with a PetscMalloc3().
101eb9baa12SBarry Smith 
102397b6df1SKris Buschelman  */
10316ebf90aSShri Abhyankar 
10416ebf90aSShri Abhyankar #undef __FUNCT__
10516ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqaij"
106bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
107b24902e0SBarry Smith {
108185f6596SHong Zhang   const PetscInt *ai,*aj,*ajj,M=A->rmap->n;
10967877ebaSShri Abhyankar   PetscInt       nz,rnz,i,j;
110dfbe8321SBarry Smith   PetscErrorCode ierr;
111c1490034SHong Zhang   PetscInt       *row,*col;
11216ebf90aSShri Abhyankar   Mat_SeqAIJ     *aa=(Mat_SeqAIJ*)A->data;
113397b6df1SKris Buschelman 
114397b6df1SKris Buschelman   PetscFunctionBegin;
11516ebf90aSShri Abhyankar   *v=aa->a;
116bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
1172205254eSKarl Rupp     nz   = aa->nz;
1182205254eSKarl Rupp     ai   = aa->i;
1192205254eSKarl Rupp     aj   = aa->j;
12016ebf90aSShri Abhyankar     *nnz = nz;
121785e854fSJed Brown     ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr);
122185f6596SHong Zhang     col  = row + nz;
123185f6596SHong Zhang 
12416ebf90aSShri Abhyankar     nz = 0;
12516ebf90aSShri Abhyankar     for (i=0; i<M; i++) {
12616ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
12767877ebaSShri Abhyankar       ajj = aj + ai[i];
12867877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
12967877ebaSShri Abhyankar         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
13016ebf90aSShri Abhyankar       }
13116ebf90aSShri Abhyankar     }
13216ebf90aSShri Abhyankar     *r = row; *c = col;
13316ebf90aSShri Abhyankar   }
13416ebf90aSShri Abhyankar   PetscFunctionReturn(0);
13516ebf90aSShri Abhyankar }
136397b6df1SKris Buschelman 
13716ebf90aSShri Abhyankar #undef __FUNCT__
13867877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqbaij_seqaij"
139bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
14067877ebaSShri Abhyankar {
14167877ebaSShri Abhyankar   Mat_SeqBAIJ    *aa=(Mat_SeqBAIJ*)A->data;
14233d57670SJed Brown   const PetscInt *ai,*aj,*ajj,bs2 = aa->bs2;
14333d57670SJed Brown   PetscInt       bs,M,nz,idx=0,rnz,i,j,k,m;
14467877ebaSShri Abhyankar   PetscErrorCode ierr;
14567877ebaSShri Abhyankar   PetscInt       *row,*col;
14667877ebaSShri Abhyankar 
14767877ebaSShri Abhyankar   PetscFunctionBegin;
14833d57670SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
14933d57670SJed Brown   M = A->rmap->N/bs;
150cf3759fdSShri Abhyankar   *v = aa->a;
151bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
152cf3759fdSShri Abhyankar     ai   = aa->i; aj = aa->j;
15367877ebaSShri Abhyankar     nz   = bs2*aa->nz;
15467877ebaSShri Abhyankar     *nnz = nz;
155785e854fSJed Brown     ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr);
156185f6596SHong Zhang     col  = row + nz;
157185f6596SHong Zhang 
15867877ebaSShri Abhyankar     for (i=0; i<M; i++) {
15967877ebaSShri Abhyankar       ajj = aj + ai[i];
16067877ebaSShri Abhyankar       rnz = ai[i+1] - ai[i];
16167877ebaSShri Abhyankar       for (k=0; k<rnz; k++) {
16267877ebaSShri Abhyankar         for (j=0; j<bs; j++) {
16367877ebaSShri Abhyankar           for (m=0; m<bs; m++) {
16467877ebaSShri Abhyankar             row[idx]   = i*bs + m + shift;
165cf3759fdSShri Abhyankar             col[idx++] = bs*(ajj[k]) + j + shift;
16667877ebaSShri Abhyankar           }
16767877ebaSShri Abhyankar         }
16867877ebaSShri Abhyankar       }
16967877ebaSShri Abhyankar     }
170cf3759fdSShri Abhyankar     *r = row; *c = col;
17167877ebaSShri Abhyankar   }
17267877ebaSShri Abhyankar   PetscFunctionReturn(0);
17367877ebaSShri Abhyankar }
17467877ebaSShri Abhyankar 
17567877ebaSShri Abhyankar #undef __FUNCT__
17616ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqsbaij_seqsbaij"
177bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
17816ebf90aSShri Abhyankar {
17967877ebaSShri Abhyankar   const PetscInt *ai, *aj,*ajj,M=A->rmap->n;
18067877ebaSShri Abhyankar   PetscInt       nz,rnz,i,j;
18116ebf90aSShri Abhyankar   PetscErrorCode ierr;
18216ebf90aSShri Abhyankar   PetscInt       *row,*col;
18316ebf90aSShri Abhyankar   Mat_SeqSBAIJ   *aa=(Mat_SeqSBAIJ*)A->data;
18416ebf90aSShri Abhyankar 
18516ebf90aSShri Abhyankar   PetscFunctionBegin;
186882afa5aSHong Zhang   *v = aa->a;
187bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
1882205254eSKarl Rupp     nz   = aa->nz;
1892205254eSKarl Rupp     ai   = aa->i;
1902205254eSKarl Rupp     aj   = aa->j;
1912205254eSKarl Rupp     *v   = aa->a;
19216ebf90aSShri Abhyankar     *nnz = nz;
193785e854fSJed Brown     ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr);
194185f6596SHong Zhang     col  = row + nz;
195185f6596SHong Zhang 
19616ebf90aSShri Abhyankar     nz = 0;
19716ebf90aSShri Abhyankar     for (i=0; i<M; i++) {
19816ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
19967877ebaSShri Abhyankar       ajj = aj + ai[i];
20067877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
20167877ebaSShri Abhyankar         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
20216ebf90aSShri Abhyankar       }
20316ebf90aSShri Abhyankar     }
20416ebf90aSShri Abhyankar     *r = row; *c = col;
20516ebf90aSShri Abhyankar   }
20616ebf90aSShri Abhyankar   PetscFunctionReturn(0);
20716ebf90aSShri Abhyankar }
20816ebf90aSShri Abhyankar 
20916ebf90aSShri Abhyankar #undef __FUNCT__
21016ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqsbaij"
211bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
21216ebf90aSShri Abhyankar {
21367877ebaSShri Abhyankar   const PetscInt    *ai,*aj,*ajj,*adiag,M=A->rmap->n;
21467877ebaSShri Abhyankar   PetscInt          nz,rnz,i,j;
21567877ebaSShri Abhyankar   const PetscScalar *av,*v1;
21616ebf90aSShri Abhyankar   PetscScalar       *val;
21716ebf90aSShri Abhyankar   PetscErrorCode    ierr;
21816ebf90aSShri Abhyankar   PetscInt          *row,*col;
219829b1710SHong Zhang   Mat_SeqAIJ        *aa=(Mat_SeqAIJ*)A->data;
22016ebf90aSShri Abhyankar 
22116ebf90aSShri Abhyankar   PetscFunctionBegin;
22216ebf90aSShri Abhyankar   ai   =aa->i; aj=aa->j;av=aa->a;
22316ebf90aSShri Abhyankar   adiag=aa->diag;
224bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
225829b1710SHong Zhang     /* count nz in the uppper triangular part of A */
226829b1710SHong Zhang     nz = 0;
227829b1710SHong Zhang     for (i=0; i<M; i++) nz += ai[i+1] - adiag[i];
22816ebf90aSShri Abhyankar     *nnz = nz;
229829b1710SHong Zhang 
230185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
231185f6596SHong Zhang     col  = row + nz;
232185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
233185f6596SHong Zhang 
23416ebf90aSShri Abhyankar     nz = 0;
23516ebf90aSShri Abhyankar     for (i=0; i<M; i++) {
23616ebf90aSShri Abhyankar       rnz = ai[i+1] - adiag[i];
23767877ebaSShri Abhyankar       ajj = aj + adiag[i];
238cf3759fdSShri Abhyankar       v1  = av + adiag[i];
23967877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
24067877ebaSShri Abhyankar         row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
24116ebf90aSShri Abhyankar       }
24216ebf90aSShri Abhyankar     }
24316ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
244397b6df1SKris Buschelman   } else {
24516ebf90aSShri Abhyankar     nz = 0; val = *v;
24616ebf90aSShri Abhyankar     for (i=0; i <M; i++) {
24716ebf90aSShri Abhyankar       rnz = ai[i+1] - adiag[i];
24867877ebaSShri Abhyankar       ajj = aj + adiag[i];
24967877ebaSShri Abhyankar       v1  = av + adiag[i];
25067877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
25167877ebaSShri Abhyankar         val[nz++] = v1[j];
25216ebf90aSShri Abhyankar       }
25316ebf90aSShri Abhyankar     }
25416ebf90aSShri Abhyankar   }
25516ebf90aSShri Abhyankar   PetscFunctionReturn(0);
25616ebf90aSShri Abhyankar }
25716ebf90aSShri Abhyankar 
25816ebf90aSShri Abhyankar #undef __FUNCT__
25916ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpisbaij_mpisbaij"
260bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
26116ebf90aSShri Abhyankar {
26216ebf90aSShri Abhyankar   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
26316ebf90aSShri Abhyankar   PetscErrorCode    ierr;
26416ebf90aSShri Abhyankar   PetscInt          rstart,nz,i,j,jj,irow,countA,countB;
26516ebf90aSShri Abhyankar   PetscInt          *row,*col;
26616ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
26716ebf90aSShri Abhyankar   PetscScalar       *val;
268397b6df1SKris Buschelman   Mat_MPISBAIJ      *mat = (Mat_MPISBAIJ*)A->data;
269397b6df1SKris Buschelman   Mat_SeqSBAIJ      *aa  = (Mat_SeqSBAIJ*)(mat->A)->data;
270397b6df1SKris Buschelman   Mat_SeqBAIJ       *bb  = (Mat_SeqBAIJ*)(mat->B)->data;
27116ebf90aSShri Abhyankar 
27216ebf90aSShri Abhyankar   PetscFunctionBegin;
273d0f46423SBarry Smith   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
274397b6df1SKris Buschelman   av=aa->a; bv=bb->a;
275397b6df1SKris Buschelman 
2762205254eSKarl Rupp   garray = mat->garray;
2772205254eSKarl Rupp 
278bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
27916ebf90aSShri Abhyankar     nz   = aa->nz + bb->nz;
28016ebf90aSShri Abhyankar     *nnz = nz;
281185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
282185f6596SHong Zhang     col  = row + nz;
283185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
284185f6596SHong Zhang 
285397b6df1SKris Buschelman     *r = row; *c = col; *v = val;
286397b6df1SKris Buschelman   } else {
287397b6df1SKris Buschelman     row = *r; col = *c; val = *v;
288397b6df1SKris Buschelman   }
289397b6df1SKris Buschelman 
290028e57e8SHong Zhang   jj = 0; irow = rstart;
291397b6df1SKris Buschelman   for (i=0; i<m; i++) {
292397b6df1SKris Buschelman     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
293397b6df1SKris Buschelman     countA = ai[i+1] - ai[i];
294397b6df1SKris Buschelman     countB = bi[i+1] - bi[i];
295397b6df1SKris Buschelman     bjj    = bj + bi[i];
29616ebf90aSShri Abhyankar     v1     = av + ai[i];
29716ebf90aSShri Abhyankar     v2     = bv + bi[i];
298397b6df1SKris Buschelman 
299397b6df1SKris Buschelman     /* A-part */
300397b6df1SKris Buschelman     for (j=0; j<countA; j++) {
301bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
302397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
303397b6df1SKris Buschelman       }
30416ebf90aSShri Abhyankar       val[jj++] = v1[j];
305397b6df1SKris Buschelman     }
30616ebf90aSShri Abhyankar 
30716ebf90aSShri Abhyankar     /* B-part */
30816ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
309bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
310397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
311397b6df1SKris Buschelman       }
31216ebf90aSShri Abhyankar       val[jj++] = v2[j];
31316ebf90aSShri Abhyankar     }
31416ebf90aSShri Abhyankar     irow++;
31516ebf90aSShri Abhyankar   }
31616ebf90aSShri Abhyankar   PetscFunctionReturn(0);
31716ebf90aSShri Abhyankar }
31816ebf90aSShri Abhyankar 
31916ebf90aSShri Abhyankar #undef __FUNCT__
32016ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpiaij"
321bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
32216ebf90aSShri Abhyankar {
32316ebf90aSShri Abhyankar   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
32416ebf90aSShri Abhyankar   PetscErrorCode    ierr;
32516ebf90aSShri Abhyankar   PetscInt          rstart,nz,i,j,jj,irow,countA,countB;
32616ebf90aSShri Abhyankar   PetscInt          *row,*col;
32716ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
32816ebf90aSShri Abhyankar   PetscScalar       *val;
32916ebf90aSShri Abhyankar   Mat_MPIAIJ        *mat = (Mat_MPIAIJ*)A->data;
33016ebf90aSShri Abhyankar   Mat_SeqAIJ        *aa  = (Mat_SeqAIJ*)(mat->A)->data;
33116ebf90aSShri Abhyankar   Mat_SeqAIJ        *bb  = (Mat_SeqAIJ*)(mat->B)->data;
33216ebf90aSShri Abhyankar 
33316ebf90aSShri Abhyankar   PetscFunctionBegin;
33416ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
33516ebf90aSShri Abhyankar   av=aa->a; bv=bb->a;
33616ebf90aSShri Abhyankar 
3372205254eSKarl Rupp   garray = mat->garray;
3382205254eSKarl Rupp 
339bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
34016ebf90aSShri Abhyankar     nz   = aa->nz + bb->nz;
34116ebf90aSShri Abhyankar     *nnz = nz;
342185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
343185f6596SHong Zhang     col  = row + nz;
344185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
345185f6596SHong Zhang 
34616ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
34716ebf90aSShri Abhyankar   } else {
34816ebf90aSShri Abhyankar     row = *r; col = *c; val = *v;
34916ebf90aSShri Abhyankar   }
35016ebf90aSShri Abhyankar 
35116ebf90aSShri Abhyankar   jj = 0; irow = rstart;
35216ebf90aSShri Abhyankar   for (i=0; i<m; i++) {
35316ebf90aSShri Abhyankar     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
35416ebf90aSShri Abhyankar     countA = ai[i+1] - ai[i];
35516ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
35616ebf90aSShri Abhyankar     bjj    = bj + bi[i];
35716ebf90aSShri Abhyankar     v1     = av + ai[i];
35816ebf90aSShri Abhyankar     v2     = bv + bi[i];
35916ebf90aSShri Abhyankar 
36016ebf90aSShri Abhyankar     /* A-part */
36116ebf90aSShri Abhyankar     for (j=0; j<countA; j++) {
362bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
36316ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
36416ebf90aSShri Abhyankar       }
36516ebf90aSShri Abhyankar       val[jj++] = v1[j];
36616ebf90aSShri Abhyankar     }
36716ebf90aSShri Abhyankar 
36816ebf90aSShri Abhyankar     /* B-part */
36916ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
370bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
37116ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
37216ebf90aSShri Abhyankar       }
37316ebf90aSShri Abhyankar       val[jj++] = v2[j];
37416ebf90aSShri Abhyankar     }
37516ebf90aSShri Abhyankar     irow++;
37616ebf90aSShri Abhyankar   }
37716ebf90aSShri Abhyankar   PetscFunctionReturn(0);
37816ebf90aSShri Abhyankar }
37916ebf90aSShri Abhyankar 
38016ebf90aSShri Abhyankar #undef __FUNCT__
38167877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpibaij_mpiaij"
382bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
38367877ebaSShri Abhyankar {
38467877ebaSShri Abhyankar   Mat_MPIBAIJ       *mat    = (Mat_MPIBAIJ*)A->data;
38567877ebaSShri Abhyankar   Mat_SeqBAIJ       *aa     = (Mat_SeqBAIJ*)(mat->A)->data;
38667877ebaSShri Abhyankar   Mat_SeqBAIJ       *bb     = (Mat_SeqBAIJ*)(mat->B)->data;
38767877ebaSShri Abhyankar   const PetscInt    *ai     = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj;
388d985c460SShri Abhyankar   const PetscInt    *garray = mat->garray,mbs=mat->mbs,rstart=A->rmap->rstart;
38933d57670SJed Brown   const PetscInt    bs2=mat->bs2;
39067877ebaSShri Abhyankar   PetscErrorCode    ierr;
39133d57670SJed Brown   PetscInt          bs,nz,i,j,k,n,jj,irow,countA,countB,idx;
39267877ebaSShri Abhyankar   PetscInt          *row,*col;
39367877ebaSShri Abhyankar   const PetscScalar *av=aa->a, *bv=bb->a,*v1,*v2;
39467877ebaSShri Abhyankar   PetscScalar       *val;
39567877ebaSShri Abhyankar 
39667877ebaSShri Abhyankar   PetscFunctionBegin;
39733d57670SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
398bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
39967877ebaSShri Abhyankar     nz   = bs2*(aa->nz + bb->nz);
40067877ebaSShri Abhyankar     *nnz = nz;
401185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
402185f6596SHong Zhang     col  = row + nz;
403185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
404185f6596SHong Zhang 
40567877ebaSShri Abhyankar     *r = row; *c = col; *v = val;
40667877ebaSShri Abhyankar   } else {
40767877ebaSShri Abhyankar     row = *r; col = *c; val = *v;
40867877ebaSShri Abhyankar   }
40967877ebaSShri Abhyankar 
410d985c460SShri Abhyankar   jj = 0; irow = rstart;
41167877ebaSShri Abhyankar   for (i=0; i<mbs; i++) {
41267877ebaSShri Abhyankar     countA = ai[i+1] - ai[i];
41367877ebaSShri Abhyankar     countB = bi[i+1] - bi[i];
41467877ebaSShri Abhyankar     ajj    = aj + ai[i];
41567877ebaSShri Abhyankar     bjj    = bj + bi[i];
41667877ebaSShri Abhyankar     v1     = av + bs2*ai[i];
41767877ebaSShri Abhyankar     v2     = bv + bs2*bi[i];
41867877ebaSShri Abhyankar 
41967877ebaSShri Abhyankar     idx = 0;
42067877ebaSShri Abhyankar     /* A-part */
42167877ebaSShri Abhyankar     for (k=0; k<countA; k++) {
42267877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
42367877ebaSShri Abhyankar         for (n=0; n<bs; n++) {
424bccb9932SShri Abhyankar           if (reuse == MAT_INITIAL_MATRIX) {
425d985c460SShri Abhyankar             row[jj] = irow + n + shift;
426d985c460SShri Abhyankar             col[jj] = rstart + bs*ajj[k] + j + shift;
42767877ebaSShri Abhyankar           }
42867877ebaSShri Abhyankar           val[jj++] = v1[idx++];
42967877ebaSShri Abhyankar         }
43067877ebaSShri Abhyankar       }
43167877ebaSShri Abhyankar     }
43267877ebaSShri Abhyankar 
43367877ebaSShri Abhyankar     idx = 0;
43467877ebaSShri Abhyankar     /* B-part */
43567877ebaSShri Abhyankar     for (k=0; k<countB; k++) {
43667877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
43767877ebaSShri Abhyankar         for (n=0; n<bs; n++) {
438bccb9932SShri Abhyankar           if (reuse == MAT_INITIAL_MATRIX) {
439d985c460SShri Abhyankar             row[jj] = irow + n + shift;
440d985c460SShri Abhyankar             col[jj] = bs*garray[bjj[k]] + j + shift;
44167877ebaSShri Abhyankar           }
442d985c460SShri Abhyankar           val[jj++] = v2[idx++];
44367877ebaSShri Abhyankar         }
44467877ebaSShri Abhyankar       }
44567877ebaSShri Abhyankar     }
446d985c460SShri Abhyankar     irow += bs;
44767877ebaSShri Abhyankar   }
44867877ebaSShri Abhyankar   PetscFunctionReturn(0);
44967877ebaSShri Abhyankar }
45067877ebaSShri Abhyankar 
45167877ebaSShri Abhyankar #undef __FUNCT__
45216ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpisbaij"
453bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
45416ebf90aSShri Abhyankar {
45516ebf90aSShri Abhyankar   const PetscInt    *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
45616ebf90aSShri Abhyankar   PetscErrorCode    ierr;
457e0bace9bSHong Zhang   PetscInt          rstart,nz,nza,nzb,i,j,jj,irow,countA,countB;
45816ebf90aSShri Abhyankar   PetscInt          *row,*col;
45916ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
46016ebf90aSShri Abhyankar   PetscScalar       *val;
46116ebf90aSShri Abhyankar   Mat_MPIAIJ        *mat =  (Mat_MPIAIJ*)A->data;
46216ebf90aSShri Abhyankar   Mat_SeqAIJ        *aa  =(Mat_SeqAIJ*)(mat->A)->data;
46316ebf90aSShri Abhyankar   Mat_SeqAIJ        *bb  =(Mat_SeqAIJ*)(mat->B)->data;
46416ebf90aSShri Abhyankar 
46516ebf90aSShri Abhyankar   PetscFunctionBegin;
46616ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j; adiag=aa->diag;
46716ebf90aSShri Abhyankar   bi=bb->i; bj=bb->j; garray = mat->garray;
46816ebf90aSShri Abhyankar   av=aa->a; bv=bb->a;
4692205254eSKarl Rupp 
47016ebf90aSShri Abhyankar   rstart = A->rmap->rstart;
47116ebf90aSShri Abhyankar 
472bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
473e0bace9bSHong Zhang     nza = 0;    /* num of upper triangular entries in mat->A, including diagonals */
474e0bace9bSHong Zhang     nzb = 0;    /* num of upper triangular entries in mat->B */
47516ebf90aSShri Abhyankar     for (i=0; i<m; i++) {
476e0bace9bSHong Zhang       nza   += (ai[i+1] - adiag[i]);
47716ebf90aSShri Abhyankar       countB = bi[i+1] - bi[i];
47816ebf90aSShri Abhyankar       bjj    = bj + bi[i];
479e0bace9bSHong Zhang       for (j=0; j<countB; j++) {
480e0bace9bSHong Zhang         if (garray[bjj[j]] > rstart) nzb++;
481e0bace9bSHong Zhang       }
482e0bace9bSHong Zhang     }
48316ebf90aSShri Abhyankar 
484e0bace9bSHong Zhang     nz   = nza + nzb; /* total nz of upper triangular part of mat */
48516ebf90aSShri Abhyankar     *nnz = nz;
486185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
487185f6596SHong Zhang     col  = row + nz;
488185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
489185f6596SHong Zhang 
49016ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
49116ebf90aSShri Abhyankar   } else {
49216ebf90aSShri Abhyankar     row = *r; col = *c; val = *v;
49316ebf90aSShri Abhyankar   }
49416ebf90aSShri Abhyankar 
49516ebf90aSShri Abhyankar   jj = 0; irow = rstart;
49616ebf90aSShri Abhyankar   for (i=0; i<m; i++) {
49716ebf90aSShri Abhyankar     ajj    = aj + adiag[i];                 /* ptr to the beginning of the diagonal of this row */
49816ebf90aSShri Abhyankar     v1     = av + adiag[i];
49916ebf90aSShri Abhyankar     countA = ai[i+1] - adiag[i];
50016ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
50116ebf90aSShri Abhyankar     bjj    = bj + bi[i];
50216ebf90aSShri Abhyankar     v2     = bv + bi[i];
50316ebf90aSShri Abhyankar 
50416ebf90aSShri Abhyankar     /* A-part */
50516ebf90aSShri Abhyankar     for (j=0; j<countA; j++) {
506bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
50716ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
50816ebf90aSShri Abhyankar       }
50916ebf90aSShri Abhyankar       val[jj++] = v1[j];
51016ebf90aSShri Abhyankar     }
51116ebf90aSShri Abhyankar 
51216ebf90aSShri Abhyankar     /* B-part */
51316ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
51416ebf90aSShri Abhyankar       if (garray[bjj[j]] > rstart) {
515bccb9932SShri Abhyankar         if (reuse == MAT_INITIAL_MATRIX) {
51616ebf90aSShri Abhyankar           row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
51716ebf90aSShri Abhyankar         }
51816ebf90aSShri Abhyankar         val[jj++] = v2[j];
51916ebf90aSShri Abhyankar       }
520397b6df1SKris Buschelman     }
521397b6df1SKris Buschelman     irow++;
522397b6df1SKris Buschelman   }
523397b6df1SKris Buschelman   PetscFunctionReturn(0);
524397b6df1SKris Buschelman }
525397b6df1SKris Buschelman 
526397b6df1SKris Buschelman #undef __FUNCT__
52720be8e61SHong Zhang #define __FUNCT__ "MatGetDiagonal_MUMPS"
52820be8e61SHong Zhang PetscErrorCode MatGetDiagonal_MUMPS(Mat A,Vec v)
52920be8e61SHong Zhang {
53020be8e61SHong Zhang   PetscFunctionBegin;
53120be8e61SHong Zhang   SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type: MUMPS factor");
53220be8e61SHong Zhang   PetscFunctionReturn(0);
53320be8e61SHong Zhang }
53420be8e61SHong Zhang 
53520be8e61SHong Zhang #undef __FUNCT__
5363924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS"
537dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A)
538dfbe8321SBarry Smith {
539a5e57a09SHong Zhang   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
540dfbe8321SBarry Smith   PetscErrorCode ierr;
541b24902e0SBarry Smith 
542397b6df1SKris Buschelman   PetscFunctionBegin;
543a5e57a09SHong Zhang   if (mumps->CleanUpMUMPS) {
544397b6df1SKris Buschelman     /* Terminate instance, deallocate memories */
545a5e57a09SHong Zhang     ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr);
546a5e57a09SHong Zhang     ierr = VecScatterDestroy(&mumps->scat_rhs);CHKERRQ(ierr);
547a5e57a09SHong Zhang     ierr = VecDestroy(&mumps->b_seq);CHKERRQ(ierr);
548*2cd7d884SHong Zhang     ierr = VecDestroy(&mumps->bb_seq);CHKERRQ(ierr);
549*2cd7d884SHong Zhang     ierr = VecDestroy(&mumps->bb_mpi);CHKERRQ(ierr);
550a5e57a09SHong Zhang     ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
551a5e57a09SHong Zhang     ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr);
552a5e57a09SHong Zhang     ierr = PetscFree(mumps->id.perm_in);CHKERRQ(ierr);
553a5e57a09SHong Zhang     ierr = PetscFree(mumps->irn);CHKERRQ(ierr);
5542205254eSKarl Rupp 
555a5e57a09SHong Zhang     mumps->id.job = JOB_END;
556a5e57a09SHong Zhang     PetscMUMPS_c(&mumps->id);
557a5e57a09SHong Zhang     ierr = MPI_Comm_free(&(mumps->comm_mumps));CHKERRQ(ierr);
558397b6df1SKris Buschelman   }
559a5e57a09SHong Zhang   if (mumps->Destroy) {
560a5e57a09SHong Zhang     ierr = (mumps->Destroy)(A);CHKERRQ(ierr);
561bf0cc555SLisandro Dalcin   }
562bf0cc555SLisandro Dalcin   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
563bf0cc555SLisandro Dalcin 
56497969023SHong Zhang   /* clear composed functions */
565bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr);
566bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetIcntl_C",NULL);CHKERRQ(ierr);
567bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetIcntl_C",NULL);CHKERRQ(ierr);
568bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetCntl_C",NULL);CHKERRQ(ierr);
569bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetCntl_C",NULL);CHKERRQ(ierr);
570bc6112feSHong Zhang 
571ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfo_C",NULL);CHKERRQ(ierr);
572ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfog_C",NULL);CHKERRQ(ierr);
573ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfo_C",NULL);CHKERRQ(ierr);
574ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfog_C",NULL);CHKERRQ(ierr);
575397b6df1SKris Buschelman   PetscFunctionReturn(0);
576397b6df1SKris Buschelman }
577397b6df1SKris Buschelman 
578397b6df1SKris Buschelman #undef __FUNCT__
579f6c57405SHong Zhang #define __FUNCT__ "MatSolve_MUMPS"
580b24902e0SBarry Smith PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x)
581b24902e0SBarry Smith {
582a5e57a09SHong Zhang   Mat_MUMPS        *mumps=(Mat_MUMPS*)A->spptr;
583d54de34fSKris Buschelman   PetscScalar      *array;
58467877ebaSShri Abhyankar   Vec              b_seq;
585329ec9b3SHong Zhang   IS               is_iden,is_petsc;
586dfbe8321SBarry Smith   PetscErrorCode   ierr;
587329ec9b3SHong Zhang   PetscInt         i;
588883f2eb9SBarry Smith   static PetscBool cite1 = PETSC_FALSE,cite2 = PETSC_FALSE;
589397b6df1SKris Buschelman 
590397b6df1SKris Buschelman   PetscFunctionBegin;
591883f2eb9SBarry 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);
592883f2eb9SBarry 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);
593a5e57a09SHong Zhang   mumps->id.nrhs = 1;
594a5e57a09SHong Zhang   b_seq          = mumps->b_seq;
595a5e57a09SHong Zhang   if (mumps->size > 1) {
596329ec9b3SHong Zhang     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
597a5e57a09SHong Zhang     ierr = VecScatterBegin(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
598a5e57a09SHong Zhang     ierr = VecScatterEnd(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
599a5e57a09SHong Zhang     if (!mumps->myid) {ierr = VecGetArray(b_seq,&array);CHKERRQ(ierr);}
600397b6df1SKris Buschelman   } else {  /* size == 1 */
601397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
602397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
603397b6df1SKris Buschelman   }
604a5e57a09SHong Zhang   if (!mumps->myid) { /* define rhs on the host */
605a5e57a09SHong Zhang     mumps->id.nrhs = 1;
606397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
6072907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
608a5e57a09SHong Zhang     mumps->id.rhs = (mumps_complex*)array;
6092907cef9SHong Zhang #else
610a5e57a09SHong Zhang     mumps->id.rhs = (mumps_double_complex*)array;
6112907cef9SHong Zhang #endif
612397b6df1SKris Buschelman #else
613a5e57a09SHong Zhang     mumps->id.rhs = array;
614397b6df1SKris Buschelman #endif
615397b6df1SKris Buschelman   }
616397b6df1SKris Buschelman 
617397b6df1SKris Buschelman   /* solve phase */
618329ec9b3SHong Zhang   /*-------------*/
619a5e57a09SHong Zhang   mumps->id.job = JOB_SOLVE;
620a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
621a5e57a09SHong 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));
622397b6df1SKris Buschelman 
623a5e57a09SHong Zhang   if (mumps->size > 1) { /* convert mumps distributed solution to petsc mpi x */
624a5e57a09SHong Zhang     if (mumps->scat_sol && mumps->ICNTL9_pre != mumps->id.ICNTL(9)) {
625a5e57a09SHong Zhang       /* when id.ICNTL(9) changes, the contents of lsol_loc may change (not its size, lsol_loc), recreates scat_sol */
626a5e57a09SHong Zhang       ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
627397b6df1SKris Buschelman     }
628a5e57a09SHong Zhang     if (!mumps->scat_sol) { /* create scatter scat_sol */
629a5e57a09SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,mumps->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */
630a5e57a09SHong Zhang       for (i=0; i<mumps->id.lsol_loc; i++) {
631a5e57a09SHong Zhang         mumps->id.isol_loc[i] -= 1; /* change Fortran style to C style */
632a5e57a09SHong Zhang       }
633a5e57a09SHong Zhang       ierr = ISCreateGeneral(PETSC_COMM_SELF,mumps->id.lsol_loc,mumps->id.isol_loc,PETSC_COPY_VALUES,&is_petsc);CHKERRQ(ierr);  /* to */
634a5e57a09SHong Zhang       ierr = VecScatterCreate(mumps->x_seq,is_iden,x,is_petsc,&mumps->scat_sol);CHKERRQ(ierr);
6356bf464f9SBarry Smith       ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
6366bf464f9SBarry Smith       ierr = ISDestroy(&is_petsc);CHKERRQ(ierr);
6372205254eSKarl Rupp 
638a5e57a09SHong Zhang       mumps->ICNTL9_pre = mumps->id.ICNTL(9); /* save current value of id.ICNTL(9) */
639397b6df1SKris Buschelman     }
640a5e57a09SHong Zhang 
641a5e57a09SHong Zhang     ierr = VecScatterBegin(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
642a5e57a09SHong Zhang     ierr = VecScatterEnd(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
643329ec9b3SHong Zhang   }
644397b6df1SKris Buschelman   PetscFunctionReturn(0);
645397b6df1SKris Buschelman }
646397b6df1SKris Buschelman 
64751d5961aSHong Zhang #undef __FUNCT__
64851d5961aSHong Zhang #define __FUNCT__ "MatSolveTranspose_MUMPS"
64951d5961aSHong Zhang PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x)
65051d5961aSHong Zhang {
651a5e57a09SHong Zhang   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
65251d5961aSHong Zhang   PetscErrorCode ierr;
65351d5961aSHong Zhang 
65451d5961aSHong Zhang   PetscFunctionBegin;
655a5e57a09SHong Zhang   mumps->id.ICNTL(9) = 0;
6562205254eSKarl Rupp 
6570ad0caddSJed Brown   ierr = MatSolve_MUMPS(A,b,x);CHKERRQ(ierr);
6582205254eSKarl Rupp 
659a5e57a09SHong Zhang   mumps->id.ICNTL(9) = 1;
66051d5961aSHong Zhang   PetscFunctionReturn(0);
66151d5961aSHong Zhang }
66251d5961aSHong Zhang 
663e0b74bf9SHong Zhang #undef __FUNCT__
664e0b74bf9SHong Zhang #define __FUNCT__ "MatMatSolve_MUMPS"
665e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X)
666e0b74bf9SHong Zhang {
667bda8bf91SBarry Smith   PetscErrorCode ierr;
668bda8bf91SBarry Smith   PetscBool      flg;
6694e34a73bSHong Zhang   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
670*2cd7d884SHong Zhang   PetscInt       i,nrhs,m,M;
671*2cd7d884SHong Zhang   PetscScalar    *array,*bray;
672bda8bf91SBarry Smith 
673e0b74bf9SHong Zhang   PetscFunctionBegin;
6740298fd71SBarry Smith   ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
675ce94432eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
6760298fd71SBarry Smith   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
677ce94432eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
6784e34a73bSHong Zhang 
679*2cd7d884SHong Zhang   ierr = MatGetLocalSize(B,&m,NULL);CHKERRQ(ierr);
680*2cd7d884SHong Zhang   ierr = MatGetSize(B,&M,&nrhs);CHKERRQ(ierr);
6814e34a73bSHong Zhang 
682*2cd7d884SHong Zhang   if (mumps->size == 1) {
683*2cd7d884SHong Zhang     /* copy B to X */
684*2cd7d884SHong Zhang     ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr);
685*2cd7d884SHong Zhang     ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr);
686*2cd7d884SHong Zhang     for (i=0; i<M*nrhs; i++) array[i] = bray[i];
687*2cd7d884SHong Zhang     ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr);
688*2cd7d884SHong Zhang 
689*2cd7d884SHong Zhang     //if (!mumps->myid) { /* define rhs on the host */
690*2cd7d884SHong Zhang     mumps->id.nrhs = nrhs;
691*2cd7d884SHong Zhang     mumps->id.lrhs = M;
692*2cd7d884SHong Zhang #if defined(PETSC_USE_COMPLEX)
693*2cd7d884SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
694*2cd7d884SHong Zhang     mumps->id.rhs = (mumps_complex*)array;
695*2cd7d884SHong Zhang #else
696*2cd7d884SHong Zhang     mumps->id.rhs = (mumps_double_complex*)array;
697*2cd7d884SHong Zhang #endif
698*2cd7d884SHong Zhang #else
699*2cd7d884SHong Zhang     mumps->id.rhs = array;
700*2cd7d884SHong Zhang #endif
701*2cd7d884SHong Zhang       //}
702*2cd7d884SHong Zhang     /* solve phase */
703*2cd7d884SHong Zhang     /*-------------*/
704*2cd7d884SHong Zhang     mumps->id.job = JOB_SOLVE;
705*2cd7d884SHong Zhang     PetscMUMPS_c(&mumps->id);
706*2cd7d884SHong 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));
707*2cd7d884SHong Zhang 
708*2cd7d884SHong Zhang     ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr);
709*2cd7d884SHong Zhang     PetscFunctionReturn(0);
710*2cd7d884SHong Zhang   }
711*2cd7d884SHong Zhang 
712*2cd7d884SHong Zhang   /* copy rhs matrix B into vector bb_mpi */
713*2cd7d884SHong Zhang   if (mumps->bb_mpi) {
714*2cd7d884SHong Zhang     ierr = VecDestroy(&mumps->bb_mpi);CHKERRQ(ierr);
715*2cd7d884SHong Zhang     ierr = VecDestroy(&mumps->bb_seq);CHKERRQ(ierr);
716*2cd7d884SHong Zhang   }
717*2cd7d884SHong Zhang   ierr = MatDenseGetArray(B,&array);CHKERRQ(ierr);
718*2cd7d884SHong Zhang   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)B),1,nrhs*m,nrhs*M,(const PetscScalar*)array,&mumps->bb_mpi);CHKERRQ(ierr);
719*2cd7d884SHong Zhang   //ierr = VecView(mumps->bb_mpi,PETSC_VIEWER_STDOUT_WORLD);
720*2cd7d884SHong Zhang   ierr = MatDenseRestoreArray(B,&array);CHKERRQ(ierr);
721*2cd7d884SHong Zhang 
722*2cd7d884SHong Zhang   /* scatter bb_mpi to bb_seq because MUMPS only supports centralized rhs */
723*2cd7d884SHong Zhang   if (!mumps->myid) {
724*2cd7d884SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,nrhs*M,&mumps->bb_seq);CHKERRQ(ierr);
725*2cd7d884SHong Zhang       //ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->N,0,1,&is_iden);CHKERRQ(ierr);
726*2cd7d884SHong Zhang     } else {
727*2cd7d884SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->bb_seq);CHKERRQ(ierr);
728*2cd7d884SHong Zhang       //ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
729*2cd7d884SHong Zhang     }
730*2cd7d884SHong Zhang 
731*2cd7d884SHong Zhang   //SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_MUMPS() is not implemented yet");
732e0b74bf9SHong Zhang   PetscFunctionReturn(0);
733e0b74bf9SHong Zhang }
734e0b74bf9SHong Zhang 
735ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX)
736a58c3f20SHong Zhang /*
737a58c3f20SHong Zhang   input:
738a58c3f20SHong Zhang    F:        numeric factor
739a58c3f20SHong Zhang   output:
740a58c3f20SHong Zhang    nneg:     total number of negative pivots
741a58c3f20SHong Zhang    nzero:    0
742a58c3f20SHong Zhang    npos:     (global dimension of F) - nneg
743a58c3f20SHong Zhang */
744a58c3f20SHong Zhang 
745a58c3f20SHong Zhang #undef __FUNCT__
746a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
747dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
748a58c3f20SHong Zhang {
749a5e57a09SHong Zhang   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
750dfbe8321SBarry Smith   PetscErrorCode ierr;
751c1490034SHong Zhang   PetscMPIInt    size;
752a58c3f20SHong Zhang 
753a58c3f20SHong Zhang   PetscFunctionBegin;
754ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)F),&size);CHKERRQ(ierr);
755bcb30aebSHong 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 */
756a5e57a09SHong 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));
757ed85ac9fSHong Zhang 
758710ac8efSHong Zhang   if (nneg) *nneg = mumps->id.INFOG(12);
759ed85ac9fSHong Zhang   if (nzero || npos) {
760ed85ac9fSHong 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");
761710ac8efSHong Zhang     if (nzero) *nzero = mumps->id.INFOG(28);
762710ac8efSHong Zhang     if (npos) *npos   = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28));
763a58c3f20SHong Zhang   }
764a58c3f20SHong Zhang   PetscFunctionReturn(0);
765a58c3f20SHong Zhang }
766ace3df97SHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */
767a58c3f20SHong Zhang 
768397b6df1SKris Buschelman #undef __FUNCT__
769f6c57405SHong Zhang #define __FUNCT__ "MatFactorNumeric_MUMPS"
7700481f469SBarry Smith PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
771af281ebdSHong Zhang {
772a5e57a09SHong Zhang   Mat_MUMPS      *mumps =(Mat_MUMPS*)(F)->spptr;
7736849ba73SBarry Smith   PetscErrorCode ierr;
774e09efc27SHong Zhang   Mat            F_diag;
775ace3abfcSBarry Smith   PetscBool      isMPIAIJ;
776397b6df1SKris Buschelman 
777397b6df1SKris Buschelman   PetscFunctionBegin;
778a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
779397b6df1SKris Buschelman 
780397b6df1SKris Buschelman   /* numerical factorization phase */
781329ec9b3SHong Zhang   /*-------------------------------*/
782a5e57a09SHong Zhang   mumps->id.job = JOB_FACTNUMERIC;
7834e34a73bSHong Zhang   if (!mumps->id.ICNTL(18)) { /* A is centralized */
784a5e57a09SHong Zhang     if (!mumps->myid) {
785397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
7862907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
787a5e57a09SHong Zhang       mumps->id.a = (mumps_complex*)mumps->val;
7882907cef9SHong Zhang #else
789a5e57a09SHong Zhang       mumps->id.a = (mumps_double_complex*)mumps->val;
7902907cef9SHong Zhang #endif
791397b6df1SKris Buschelman #else
792a5e57a09SHong Zhang       mumps->id.a = mumps->val;
793397b6df1SKris Buschelman #endif
794397b6df1SKris Buschelman     }
795397b6df1SKris Buschelman   } else {
796397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
7972907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
798a5e57a09SHong Zhang     mumps->id.a_loc = (mumps_complex*)mumps->val;
7992907cef9SHong Zhang #else
800a5e57a09SHong Zhang     mumps->id.a_loc = (mumps_double_complex*)mumps->val;
8012907cef9SHong Zhang #endif
802397b6df1SKris Buschelman #else
803a5e57a09SHong Zhang     mumps->id.a_loc = mumps->val;
804397b6df1SKris Buschelman #endif
805397b6df1SKris Buschelman   }
806a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
807a5e57a09SHong Zhang   if (mumps->id.INFOG(1) < 0) {
808151787a6SHong Zhang     if (mumps->id.INFO(1) == -13) {
809151787a6SHong Zhang       if (mumps->id.INFO(2) < 0) {
810151787a6SHong 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));
811151787a6SHong Zhang       } else {
812151787a6SHong 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));
813151787a6SHong Zhang       }
814151787a6SHong 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));
815397b6df1SKris Buschelman   }
816a5e57a09SHong 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));
817397b6df1SKris Buschelman 
818dcd589f8SShri Abhyankar   (F)->assembled      = PETSC_TRUE;
819a5e57a09SHong Zhang   mumps->matstruc     = SAME_NONZERO_PATTERN;
820a5e57a09SHong Zhang   mumps->CleanUpMUMPS = PETSC_TRUE;
82167877ebaSShri Abhyankar 
822a5e57a09SHong Zhang   if (mumps->size > 1) {
82367877ebaSShri Abhyankar     PetscInt    lsol_loc;
82467877ebaSShri Abhyankar     PetscScalar *sol_loc;
8252205254eSKarl Rupp 
826c2093ab7SHong Zhang     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr);
827c2093ab7SHong Zhang     if (isMPIAIJ) F_diag = ((Mat_MPIAIJ*)(F)->data)->A;
828c2093ab7SHong Zhang     else F_diag = ((Mat_MPISBAIJ*)(F)->data)->A;
829c2093ab7SHong Zhang     F_diag->assembled = PETSC_TRUE;
830c2093ab7SHong Zhang 
831c2093ab7SHong Zhang     /* distributed solution; Create x_seq=sol_loc for repeated use */
832c2093ab7SHong Zhang     if (mumps->x_seq) {
833c2093ab7SHong Zhang       ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
834c2093ab7SHong Zhang       ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr);
835c2093ab7SHong Zhang       ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr);
836c2093ab7SHong Zhang     }
837a5e57a09SHong Zhang     lsol_loc = mumps->id.INFO(23); /* length of sol_loc */
838dcca6d9dSJed Brown     ierr = PetscMalloc2(lsol_loc,&sol_loc,lsol_loc,&mumps->id.isol_loc);CHKERRQ(ierr);
839a5e57a09SHong Zhang     mumps->id.lsol_loc = lsol_loc;
84067877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
8412907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
842a5e57a09SHong Zhang     mumps->id.sol_loc = (mumps_complex*)sol_loc;
8432907cef9SHong Zhang #else
844a5e57a09SHong Zhang     mumps->id.sol_loc = (mumps_double_complex*)sol_loc;
8452907cef9SHong Zhang #endif
84667877ebaSShri Abhyankar #else
847a5e57a09SHong Zhang     mumps->id.sol_loc = sol_loc;
84867877ebaSShri Abhyankar #endif
849a5e57a09SHong Zhang     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq);CHKERRQ(ierr);
85067877ebaSShri Abhyankar   }
851397b6df1SKris Buschelman   PetscFunctionReturn(0);
852397b6df1SKris Buschelman }
853397b6df1SKris Buschelman 
8549a2535b5SHong Zhang /* Sets MUMPS options from the options database */
855dcd589f8SShri Abhyankar #undef __FUNCT__
8569a2535b5SHong Zhang #define __FUNCT__ "PetscSetMUMPSFromOptions"
8579a2535b5SHong Zhang PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A)
858dcd589f8SShri Abhyankar {
8599a2535b5SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
860dcd589f8SShri Abhyankar   PetscErrorCode ierr;
861dcd589f8SShri Abhyankar   PetscInt       icntl;
862ace3abfcSBarry Smith   PetscBool      flg;
863dcd589f8SShri Abhyankar 
864dcd589f8SShri Abhyankar   PetscFunctionBegin;
865ce94432eSBarry Smith   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
8669a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);CHKERRQ(ierr);
8679a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(1) = icntl;
8689a2535b5SHong 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);
8699a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(2) = icntl;
8709a2535b5SHong 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);
8719a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(3) = icntl;
872dcd589f8SShri Abhyankar 
8739a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",mumps->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
8749a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(4) = icntl;
8759a2535b5SHong Zhang   if (mumps->id.ICNTL(4) || PetscLogPrintInfo) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */
8769a2535b5SHong Zhang 
877d341cd04SHong 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);
8789a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(6) = icntl;
8799a2535b5SHong Zhang 
880d341cd04SHong 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);
881dcd589f8SShri Abhyankar   if (flg) {
8822205254eSKarl 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");
8832205254eSKarl Rupp     else mumps->id.ICNTL(7) = icntl;
884dcd589f8SShri Abhyankar   }
885e0b74bf9SHong Zhang 
8860298fd71SBarry 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);
887d341cd04SHong 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() */
8880298fd71SBarry 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);
889d341cd04SHong 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);
890d341cd04SHong 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);
891d341cd04SHong 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);
892d341cd04SHong 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);
893d341cd04SHong 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);
8944e34a73bSHong 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 */
895d341cd04SHong 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 */
8969a2535b5SHong Zhang 
897d341cd04SHong 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);
8980298fd71SBarry 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);
8990298fd71SBarry 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);
9009a2535b5SHong Zhang   if (mumps->id.ICNTL(24)) {
9019a2535b5SHong Zhang     mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */
902d7ebd59bSHong Zhang   }
903d7ebd59bSHong Zhang 
904d341cd04SHong 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);
905d341cd04SHong 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);
906*2cd7d884SHong 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);
9070298fd71SBarry 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);
908d341cd04SHong 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);
9090298fd71SBarry 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);
910d341cd04SHong 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);
9114e34a73bSHong 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 */
9120298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),NULL);CHKERRQ(ierr);
913dcd589f8SShri Abhyankar 
9140298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),NULL);CHKERRQ(ierr);
9150298fd71SBarry 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);
9160298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),NULL);CHKERRQ(ierr);
9170298fd71SBarry 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);
9180298fd71SBarry 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);
919e5bb22a1SHong Zhang 
9200298fd71SBarry Smith   ierr = PetscOptionsString("-mat_mumps_ooc_tmpdir", "out of core directory", "None", mumps->id.ooc_tmpdir, mumps->id.ooc_tmpdir, 256, NULL);
921dcd589f8SShri Abhyankar   PetscOptionsEnd();
922dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
923dcd589f8SShri Abhyankar }
924dcd589f8SShri Abhyankar 
925dcd589f8SShri Abhyankar #undef __FUNCT__
926dcd589f8SShri Abhyankar #define __FUNCT__ "PetscInitializeMUMPS"
927f697e70eSHong Zhang PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps)
928dcd589f8SShri Abhyankar {
929dcd589f8SShri Abhyankar   PetscErrorCode ierr;
930dcd589f8SShri Abhyankar 
931dcd589f8SShri Abhyankar   PetscFunctionBegin;
932ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)A), &mumps->myid);
933ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&mumps->size);CHKERRQ(ierr);
934ce94432eSBarry Smith   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(mumps->comm_mumps));CHKERRQ(ierr);
9352205254eSKarl Rupp 
936f697e70eSHong Zhang   mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps);
937f697e70eSHong Zhang 
938f697e70eSHong Zhang   mumps->id.job = JOB_INIT;
939f697e70eSHong Zhang   mumps->id.par = 1;  /* host participates factorizaton and solve */
940f697e70eSHong Zhang   mumps->id.sym = mumps->sym;
9412907cef9SHong Zhang   PetscMUMPS_c(&mumps->id);
942f697e70eSHong Zhang 
943f697e70eSHong Zhang   mumps->CleanUpMUMPS = PETSC_FALSE;
9440298fd71SBarry Smith   mumps->scat_rhs     = NULL;
9450298fd71SBarry Smith   mumps->scat_sol     = NULL;
9469a2535b5SHong Zhang 
94770544d5fSHong Zhang   /* set PETSc-MUMPS default options - override MUMPS default */
9489a2535b5SHong Zhang   mumps->id.ICNTL(3) = 0;
9499a2535b5SHong Zhang   mumps->id.ICNTL(4) = 0;
9509a2535b5SHong Zhang   if (mumps->size == 1) {
9519a2535b5SHong Zhang     mumps->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
9529a2535b5SHong Zhang   } else {
9539a2535b5SHong Zhang     mumps->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
9544e34a73bSHong Zhang     mumps->id.ICNTL(20) = 0;   /* rhs is in dense format */
95570544d5fSHong Zhang     mumps->id.ICNTL(21) = 1;   /* distributed solution */
9569a2535b5SHong Zhang   }
957dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
958dcd589f8SShri Abhyankar }
959dcd589f8SShri Abhyankar 
960a5e57a09SHong 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 */
961397b6df1SKris Buschelman #undef __FUNCT__
962f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
9630481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
964b24902e0SBarry Smith {
965a5e57a09SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
966dcd589f8SShri Abhyankar   PetscErrorCode ierr;
96767877ebaSShri Abhyankar   Vec            b;
96867877ebaSShri Abhyankar   IS             is_iden;
96967877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
970397b6df1SKris Buschelman 
971397b6df1SKris Buschelman   PetscFunctionBegin;
972a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
973dcd589f8SShri Abhyankar 
9749a2535b5SHong Zhang   /* Set MUMPS options from the options database */
9759a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
976dcd589f8SShri Abhyankar 
977a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
978dcd589f8SShri Abhyankar 
97967877ebaSShri Abhyankar   /* analysis phase */
98067877ebaSShri Abhyankar   /*----------------*/
981a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
982a5e57a09SHong Zhang   mumps->id.n   = M;
983a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
98467877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
985a5e57a09SHong Zhang     if (!mumps->myid) {
986a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
987a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
98867877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
9892907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
990a5e57a09SHong Zhang         mumps->id.a = (mumps_complex*)mumps->val;
9912907cef9SHong Zhang #else
992a5e57a09SHong Zhang         mumps->id.a = (mumps_double_complex*)mumps->val;
9932907cef9SHong Zhang #endif
99467877ebaSShri Abhyankar #else
995a5e57a09SHong Zhang         mumps->id.a = mumps->val;
99667877ebaSShri Abhyankar #endif
99767877ebaSShri Abhyankar       }
998a5e57a09SHong Zhang       if (mumps->id.ICNTL(7) == 1) { /* use user-provide matrix ordering - assuming r = c ordering */
9995248a706SHong Zhang         /*
10005248a706SHong Zhang         PetscBool      flag;
10015248a706SHong Zhang         ierr = ISEqual(r,c,&flag);CHKERRQ(ierr);
10025248a706SHong Zhang         if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"row_perm != col_perm");
10035248a706SHong Zhang         ierr = ISView(r,PETSC_VIEWER_STDOUT_SELF);
10045248a706SHong Zhang          */
1005a5e57a09SHong Zhang         if (!mumps->myid) {
1006e0b74bf9SHong Zhang           const PetscInt *idx;
1007e0b74bf9SHong Zhang           PetscInt       i,*perm_in;
10082205254eSKarl Rupp 
1009785e854fSJed Brown           ierr = PetscMalloc1(M,&perm_in);CHKERRQ(ierr);
1010e0b74bf9SHong Zhang           ierr = ISGetIndices(r,&idx);CHKERRQ(ierr);
10112205254eSKarl Rupp 
1012a5e57a09SHong Zhang           mumps->id.perm_in = perm_in;
1013e0b74bf9SHong Zhang           for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */
1014e0b74bf9SHong Zhang           ierr = ISRestoreIndices(r,&idx);CHKERRQ(ierr);
1015e0b74bf9SHong Zhang         }
1016e0b74bf9SHong Zhang       }
101767877ebaSShri Abhyankar     }
101867877ebaSShri Abhyankar     break;
101967877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1020a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1021a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1022a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
102367877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
10242907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1025a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_complex*)mumps->val;
10262907cef9SHong Zhang #else
1027a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_double_complex*)mumps->val;
10282907cef9SHong Zhang #endif
102967877ebaSShri Abhyankar #else
1030a5e57a09SHong Zhang       mumps->id.a_loc = mumps->val;
103167877ebaSShri Abhyankar #endif
103267877ebaSShri Abhyankar     }
103367877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1034a5e57a09SHong Zhang     if (!mumps->myid) {
1035*2cd7d884SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->rmap->N,&mumps->b_seq);CHKERRQ(ierr);
1036*2cd7d884SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->N,0,1,&is_iden);CHKERRQ(ierr);
103767877ebaSShri Abhyankar     } else {
1038a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
103967877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
104067877ebaSShri Abhyankar     }
10412a7a6963SBarry Smith     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1042a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
10436bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
10446bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
104567877ebaSShri Abhyankar     break;
104667877ebaSShri Abhyankar   }
1047a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
1048a5e57a09SHong 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));
104967877ebaSShri Abhyankar 
1050719d5645SBarry Smith   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1051dcd589f8SShri Abhyankar   F->ops->solve           = MatSolve_MUMPS;
105251d5961aSHong Zhang   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
10534e34a73bSHong Zhang   F->ops->matsolve        = MatMatSolve_MUMPS;
1054b24902e0SBarry Smith   PetscFunctionReturn(0);
1055b24902e0SBarry Smith }
1056b24902e0SBarry Smith 
1057450b117fSShri Abhyankar /* Note the Petsc r and c permutations are ignored */
1058450b117fSShri Abhyankar #undef __FUNCT__
1059450b117fSShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_BAIJMUMPS"
1060450b117fSShri Abhyankar PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1061450b117fSShri Abhyankar {
1062a5e57a09SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
1063dcd589f8SShri Abhyankar   PetscErrorCode ierr;
106467877ebaSShri Abhyankar   Vec            b;
106567877ebaSShri Abhyankar   IS             is_iden;
106667877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1067450b117fSShri Abhyankar 
1068450b117fSShri Abhyankar   PetscFunctionBegin;
1069a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1070dcd589f8SShri Abhyankar 
10719a2535b5SHong Zhang   /* Set MUMPS options from the options database */
10729a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1073dcd589f8SShri Abhyankar 
1074a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
107567877ebaSShri Abhyankar 
107667877ebaSShri Abhyankar   /* analysis phase */
107767877ebaSShri Abhyankar   /*----------------*/
1078a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1079a5e57a09SHong Zhang   mumps->id.n   = M;
1080a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
108167877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1082a5e57a09SHong Zhang     if (!mumps->myid) {
1083a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1084a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
108567877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
10862907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1087a5e57a09SHong Zhang         mumps->id.a = (mumps_complex*)mumps->val;
10882907cef9SHong Zhang #else
1089a5e57a09SHong Zhang         mumps->id.a = (mumps_double_complex*)mumps->val;
10902907cef9SHong Zhang #endif
109167877ebaSShri Abhyankar #else
1092a5e57a09SHong Zhang         mumps->id.a = mumps->val;
109367877ebaSShri Abhyankar #endif
109467877ebaSShri Abhyankar       }
109567877ebaSShri Abhyankar     }
109667877ebaSShri Abhyankar     break;
109767877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1098a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1099a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1100a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
110167877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
11022907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1103a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_complex*)mumps->val;
11042907cef9SHong Zhang #else
1105a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_double_complex*)mumps->val;
11062907cef9SHong Zhang #endif
110767877ebaSShri Abhyankar #else
1108a5e57a09SHong Zhang       mumps->id.a_loc = mumps->val;
110967877ebaSShri Abhyankar #endif
111067877ebaSShri Abhyankar     }
111167877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1112a5e57a09SHong Zhang     if (!mumps->myid) {
1113a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
111467877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
111567877ebaSShri Abhyankar     } else {
1116a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
111767877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
111867877ebaSShri Abhyankar     }
11192a7a6963SBarry Smith     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1120a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
11216bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
11226bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
112367877ebaSShri Abhyankar     break;
112467877ebaSShri Abhyankar   }
1125a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
1126a5e57a09SHong 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));
112767877ebaSShri Abhyankar 
1128450b117fSShri Abhyankar   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1129dcd589f8SShri Abhyankar   F->ops->solve           = MatSolve_MUMPS;
113051d5961aSHong Zhang   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1131450b117fSShri Abhyankar   PetscFunctionReturn(0);
1132450b117fSShri Abhyankar }
1133b24902e0SBarry Smith 
1134141f4205SHong Zhang /* Note the Petsc r permutation and factor info are ignored */
1135397b6df1SKris Buschelman #undef __FUNCT__
113667877ebaSShri Abhyankar #define __FUNCT__ "MatCholeskyFactorSymbolic_MUMPS"
113767877ebaSShri Abhyankar PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
1138b24902e0SBarry Smith {
1139a5e57a09SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
1140dcd589f8SShri Abhyankar   PetscErrorCode ierr;
114167877ebaSShri Abhyankar   Vec            b;
114267877ebaSShri Abhyankar   IS             is_iden;
114367877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1144397b6df1SKris Buschelman 
1145397b6df1SKris Buschelman   PetscFunctionBegin;
1146a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1147dcd589f8SShri Abhyankar 
11489a2535b5SHong Zhang   /* Set MUMPS options from the options database */
11499a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1150dcd589f8SShri Abhyankar 
1151a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1152dcd589f8SShri Abhyankar 
115367877ebaSShri Abhyankar   /* analysis phase */
115467877ebaSShri Abhyankar   /*----------------*/
1155a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1156a5e57a09SHong Zhang   mumps->id.n   = M;
1157a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
115867877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1159a5e57a09SHong Zhang     if (!mumps->myid) {
1160a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1161a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
116267877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
11632907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1164a5e57a09SHong Zhang         mumps->id.a = (mumps_complex*)mumps->val;
11652907cef9SHong Zhang #else
1166a5e57a09SHong Zhang         mumps->id.a = (mumps_double_complex*)mumps->val;
11672907cef9SHong Zhang #endif
116867877ebaSShri Abhyankar #else
1169a5e57a09SHong Zhang         mumps->id.a = mumps->val;
117067877ebaSShri Abhyankar #endif
117167877ebaSShri Abhyankar       }
117267877ebaSShri Abhyankar     }
117367877ebaSShri Abhyankar     break;
117467877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1175a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1176a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1177a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
117867877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
11792907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1180a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_complex*)mumps->val;
11812907cef9SHong Zhang #else
1182a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_double_complex*)mumps->val;
11832907cef9SHong Zhang #endif
118467877ebaSShri Abhyankar #else
1185a5e57a09SHong Zhang       mumps->id.a_loc = mumps->val;
118667877ebaSShri Abhyankar #endif
118767877ebaSShri Abhyankar     }
118867877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1189a5e57a09SHong Zhang     if (!mumps->myid) {
1190a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
119167877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
119267877ebaSShri Abhyankar     } else {
1193a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
119467877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
119567877ebaSShri Abhyankar     }
11962a7a6963SBarry Smith     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1197a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
11986bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
11996bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
120067877ebaSShri Abhyankar     break;
120167877ebaSShri Abhyankar   }
1202a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
1203a5e57a09SHong 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));
120467877ebaSShri Abhyankar 
12052792810eSHong Zhang   F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
1206dcd589f8SShri Abhyankar   F->ops->solve                 = MatSolve_MUMPS;
120751d5961aSHong Zhang   F->ops->solvetranspose        = MatSolve_MUMPS;
12084e34a73bSHong Zhang   F->ops->matsolve              = MatMatSolve_MUMPS;
12094e34a73bSHong Zhang #if defined(PETSC_USE_COMPLEX)
12100298fd71SBarry Smith   F->ops->getinertia = NULL;
12114e34a73bSHong Zhang #else
12124e34a73bSHong Zhang   F->ops->getinertia = MatGetInertia_SBAIJMUMPS;
1213db4efbfdSBarry Smith #endif
1214b24902e0SBarry Smith   PetscFunctionReturn(0);
1215b24902e0SBarry Smith }
1216b24902e0SBarry Smith 
12174e34a73bSHong Zhang //update!!!
1218397b6df1SKris Buschelman #undef __FUNCT__
121964e6c443SBarry Smith #define __FUNCT__ "MatView_MUMPS"
122064e6c443SBarry Smith PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
122174ed9c26SBarry Smith {
1222f6c57405SHong Zhang   PetscErrorCode    ierr;
122364e6c443SBarry Smith   PetscBool         iascii;
122464e6c443SBarry Smith   PetscViewerFormat format;
1225a5e57a09SHong Zhang   Mat_MUMPS         *mumps=(Mat_MUMPS*)A->spptr;
1226f6c57405SHong Zhang 
1227f6c57405SHong Zhang   PetscFunctionBegin;
122864e6c443SBarry Smith   /* check if matrix is mumps type */
122964e6c443SBarry Smith   if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0);
123064e6c443SBarry Smith 
1231251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
123264e6c443SBarry Smith   if (iascii) {
123364e6c443SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
123464e6c443SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO) {
123564e6c443SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
1236a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                   %d \n",mumps->id.sym);CHKERRQ(ierr);
1237a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):            %d \n",mumps->id.par);CHKERRQ(ierr);
1238a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):         %d \n",mumps->id.ICNTL(1));CHKERRQ(ierr);
1239a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg): %d \n",mumps->id.ICNTL(2));CHKERRQ(ierr);
1240a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):   %d \n",mumps->id.ICNTL(3));CHKERRQ(ierr);
1241a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):        %d \n",mumps->id.ICNTL(4));CHKERRQ(ierr);
1242a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):         %d \n",mumps->id.ICNTL(5));CHKERRQ(ierr);
1243a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):        %d \n",mumps->id.ICNTL(6));CHKERRQ(ierr);
1244a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (sequentia matrix ordering):%d \n",mumps->id.ICNTL(7));CHKERRQ(ierr);
1245a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scalling strategy):        %d \n",mumps->id.ICNTL(8));CHKERRQ(ierr);
1246a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements):  %d \n",mumps->id.ICNTL(10));CHKERRQ(ierr);
1247a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):          %d \n",mumps->id.ICNTL(11));CHKERRQ(ierr);
1248a5e57a09SHong Zhang       if (mumps->id.ICNTL(11)>0) {
1249a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(4) (inf norm of input mat):        %g\n",mumps->id.RINFOG(4));CHKERRQ(ierr);
1250a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(5) (inf norm of solution):         %g\n",mumps->id.RINFOG(5));CHKERRQ(ierr);
1251a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(6) (inf norm of residual):         %g\n",mumps->id.RINFOG(6));CHKERRQ(ierr);
1252a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8));CHKERRQ(ierr);
1253a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(9) (error estimate):               %g \n",mumps->id.RINFOG(9));CHKERRQ(ierr);
1254a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11));CHKERRQ(ierr);
1255f6c57405SHong Zhang       }
1256a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",mumps->id.ICNTL(12));CHKERRQ(ierr);
1257a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",mumps->id.ICNTL(13));CHKERRQ(ierr);
1258a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",mumps->id.ICNTL(14));CHKERRQ(ierr);
1259f6c57405SHong Zhang       /* ICNTL(15-17) not used */
1260a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",mumps->id.ICNTL(18));CHKERRQ(ierr);
1261a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Shur complement info):                       %d \n",mumps->id.ICNTL(19));CHKERRQ(ierr);
1262a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",mumps->id.ICNTL(20));CHKERRQ(ierr);
1263a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (somumpstion struct):                            %d \n",mumps->id.ICNTL(21));CHKERRQ(ierr);
1264a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",mumps->id.ICNTL(22));CHKERRQ(ierr);
1265a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",mumps->id.ICNTL(23));CHKERRQ(ierr);
1266c0165424SHong Zhang 
1267a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",mumps->id.ICNTL(24));CHKERRQ(ierr);
1268a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",mumps->id.ICNTL(25));CHKERRQ(ierr);
1269a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",mumps->id.ICNTL(26));CHKERRQ(ierr);
1270a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",mumps->id.ICNTL(27));CHKERRQ(ierr);
1271a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(28) (use parallel or sequential ordering):        %d \n",mumps->id.ICNTL(28));CHKERRQ(ierr);
1272a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(29) (parallel ordering):                          %d \n",mumps->id.ICNTL(29));CHKERRQ(ierr);
127342179a6aSHong Zhang 
1274a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(30) (user-specified set of entries in inv(A)):    %d \n",mumps->id.ICNTL(30));CHKERRQ(ierr);
1275a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(31) (factors is discarded in the solve phase):    %d \n",mumps->id.ICNTL(31));CHKERRQ(ierr);
1276a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(33) (compute determinant):                        %d \n",mumps->id.ICNTL(33));CHKERRQ(ierr);
1277f6c57405SHong Zhang 
1278a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",mumps->id.CNTL(1));CHKERRQ(ierr);
1279a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2));CHKERRQ(ierr);
1280a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absomumpste pivoting threshold):      %g \n",mumps->id.CNTL(3));CHKERRQ(ierr);
1281a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(4) (vamumpse of static pivoting):         %g \n",mumps->id.CNTL(4));CHKERRQ(ierr);
1282a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",mumps->id.CNTL(5));CHKERRQ(ierr);
1283f6c57405SHong Zhang 
1284f6c57405SHong Zhang       /* infomation local to each processor */
128534ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);
12867b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
1287a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %g \n",mumps->myid,mumps->id.RINFO(1));CHKERRQ(ierr);
128834ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
128934ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);
1290a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(2));CHKERRQ(ierr);
129134ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
129234ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);
1293a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(3));CHKERRQ(ierr);
129434ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1295f6c57405SHong Zhang 
129634ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr);
1297a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] %d \n",mumps->myid,mumps->id.INFO(15));CHKERRQ(ierr);
129834ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1299f6c57405SHong Zhang 
130034ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr);
1301a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(16));CHKERRQ(ierr);
130234ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1303f6c57405SHong Zhang 
130434ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);
1305a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(23));CHKERRQ(ierr);
130634ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
13077b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
1308f6c57405SHong Zhang 
1309a5e57a09SHong Zhang       if (!mumps->myid) { /* information from the host */
1310a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",mumps->id.RINFOG(1));CHKERRQ(ierr);
1311a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",mumps->id.RINFOG(2));CHKERRQ(ierr);
1312a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",mumps->id.RINFOG(3));CHKERRQ(ierr);
1313a5e57a09SHong 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);
1314f6c57405SHong Zhang 
1315a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(3));CHKERRQ(ierr);
1316a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(4));CHKERRQ(ierr);
1317a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",mumps->id.INFOG(5));CHKERRQ(ierr);
1318a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",mumps->id.INFOG(6));CHKERRQ(ierr);
1319a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively use after analysis): %d \n",mumps->id.INFOG(7));CHKERRQ(ierr);
1320a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",mumps->id.INFOG(8));CHKERRQ(ierr);
1321a5e57a09SHong 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);
1322a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",mumps->id.INFOG(10));CHKERRQ(ierr);
1323a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",mumps->id.INFOG(11));CHKERRQ(ierr);
1324a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",mumps->id.INFOG(12));CHKERRQ(ierr);
1325a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",mumps->id.INFOG(13));CHKERRQ(ierr);
1326a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",mumps->id.INFOG(14));CHKERRQ(ierr);
1327a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",mumps->id.INFOG(15));CHKERRQ(ierr);
1328a5e57a09SHong 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);
1329a5e57a09SHong 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);
1330a5e57a09SHong 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);
1331a5e57a09SHong 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);
1332a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",mumps->id.INFOG(20));CHKERRQ(ierr);
1333a5e57a09SHong 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);
1334a5e57a09SHong 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);
1335a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",mumps->id.INFOG(23));CHKERRQ(ierr);
1336a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",mumps->id.INFOG(24));CHKERRQ(ierr);
1337a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",mumps->id.INFOG(25));CHKERRQ(ierr);
133840d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(28) (after factorization: number of null pivots encountered): %d\n",mumps->id.INFOG(28));CHKERRQ(ierr);
133940d435e3SHong 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);
134040d435e3SHong 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);
134140d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(32) (after analysis: type of analysis done): %d\n",mumps->id.INFOG(32));CHKERRQ(ierr);
134240d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(33) (value used for ICNTL(8)): %d\n",mumps->id.INFOG(33));CHKERRQ(ierr);
134340d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(34) (exponent of the determinant if determinant is requested): %d\n",mumps->id.INFOG(34));CHKERRQ(ierr);
1344f6c57405SHong Zhang       }
1345f6c57405SHong Zhang     }
1346cb828f0fSHong Zhang   }
1347f6c57405SHong Zhang   PetscFunctionReturn(0);
1348f6c57405SHong Zhang }
1349f6c57405SHong Zhang 
135035bd34faSBarry Smith #undef __FUNCT__
135135bd34faSBarry Smith #define __FUNCT__ "MatGetInfo_MUMPS"
135235bd34faSBarry Smith PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
135335bd34faSBarry Smith {
1354cb828f0fSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)A->spptr;
135535bd34faSBarry Smith 
135635bd34faSBarry Smith   PetscFunctionBegin;
135735bd34faSBarry Smith   info->block_size        = 1.0;
1358cb828f0fSHong Zhang   info->nz_allocated      = mumps->id.INFOG(20);
1359cb828f0fSHong Zhang   info->nz_used           = mumps->id.INFOG(20);
136035bd34faSBarry Smith   info->nz_unneeded       = 0.0;
136135bd34faSBarry Smith   info->assemblies        = 0.0;
136235bd34faSBarry Smith   info->mallocs           = 0.0;
136335bd34faSBarry Smith   info->memory            = 0.0;
136435bd34faSBarry Smith   info->fill_ratio_given  = 0;
136535bd34faSBarry Smith   info->fill_ratio_needed = 0;
136635bd34faSBarry Smith   info->factor_mallocs    = 0;
136735bd34faSBarry Smith   PetscFunctionReturn(0);
136835bd34faSBarry Smith }
136935bd34faSBarry Smith 
13705ccb76cbSHong Zhang /* -------------------------------------------------------------------------------------------*/
13715ccb76cbSHong Zhang #undef __FUNCT__
13725ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl_MUMPS"
13735ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival)
13745ccb76cbSHong Zhang {
1375a5e57a09SHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
13765ccb76cbSHong Zhang 
13775ccb76cbSHong Zhang   PetscFunctionBegin;
1378a5e57a09SHong Zhang   mumps->id.ICNTL(icntl) = ival;
13795ccb76cbSHong Zhang   PetscFunctionReturn(0);
13805ccb76cbSHong Zhang }
13815ccb76cbSHong Zhang 
13825ccb76cbSHong Zhang #undef __FUNCT__
1383bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetIcntl_MUMPS"
1384bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt *ival)
1385bc6112feSHong Zhang {
1386bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1387bc6112feSHong Zhang 
1388bc6112feSHong Zhang   PetscFunctionBegin;
1389bc6112feSHong Zhang   *ival = mumps->id.ICNTL(icntl);
1390bc6112feSHong Zhang   PetscFunctionReturn(0);
1391bc6112feSHong Zhang }
1392bc6112feSHong Zhang 
1393bc6112feSHong Zhang #undef __FUNCT__
13945ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl"
13955ccb76cbSHong Zhang /*@
13965ccb76cbSHong Zhang   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()
13975ccb76cbSHong Zhang 
13985ccb76cbSHong Zhang    Logically Collective on Mat
13995ccb76cbSHong Zhang 
14005ccb76cbSHong Zhang    Input Parameters:
14015ccb76cbSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
14025ccb76cbSHong Zhang .  icntl - index of MUMPS parameter array ICNTL()
14035ccb76cbSHong Zhang -  ival - value of MUMPS ICNTL(icntl)
14045ccb76cbSHong Zhang 
14055ccb76cbSHong Zhang   Options Database:
14065ccb76cbSHong Zhang .   -mat_mumps_icntl_<icntl> <ival>
14075ccb76cbSHong Zhang 
14085ccb76cbSHong Zhang    Level: beginner
14095ccb76cbSHong Zhang 
14105ccb76cbSHong Zhang    References: MUMPS Users' Guide
14115ccb76cbSHong Zhang 
14125ccb76cbSHong Zhang .seealso: MatGetFactor()
14135ccb76cbSHong Zhang @*/
14145ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival)
14155ccb76cbSHong Zhang {
14165ccb76cbSHong Zhang   PetscErrorCode ierr;
14175ccb76cbSHong Zhang 
14185ccb76cbSHong Zhang   PetscFunctionBegin;
14195ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
14205ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,ival,3);
14215ccb76cbSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr);
14225ccb76cbSHong Zhang   PetscFunctionReturn(0);
14235ccb76cbSHong Zhang }
14245ccb76cbSHong Zhang 
1425bc6112feSHong Zhang #undef __FUNCT__
1426bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetIcntl"
1427a21f80fcSHong Zhang /*@
1428a21f80fcSHong Zhang   MatMumpsGetIcntl - Get MUMPS parameter ICNTL()
1429a21f80fcSHong Zhang 
1430a21f80fcSHong Zhang    Logically Collective on Mat
1431a21f80fcSHong Zhang 
1432a21f80fcSHong Zhang    Input Parameters:
1433a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1434a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array ICNTL()
1435a21f80fcSHong Zhang 
1436a21f80fcSHong Zhang   Output Parameter:
1437a21f80fcSHong Zhang .  ival - value of MUMPS ICNTL(icntl)
1438a21f80fcSHong Zhang 
1439a21f80fcSHong Zhang    Level: beginner
1440a21f80fcSHong Zhang 
1441a21f80fcSHong Zhang    References: MUMPS Users' Guide
1442a21f80fcSHong Zhang 
1443a21f80fcSHong Zhang .seealso: MatGetFactor()
1444a21f80fcSHong Zhang @*/
1445bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl(Mat F,PetscInt icntl,PetscInt *ival)
1446bc6112feSHong Zhang {
1447bc6112feSHong Zhang   PetscErrorCode ierr;
1448bc6112feSHong Zhang 
1449bc6112feSHong Zhang   PetscFunctionBegin;
1450bc6112feSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
1451bc6112feSHong Zhang   PetscValidIntPointer(ival,3);
1452bc6112feSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetIcntl_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
1453bc6112feSHong Zhang   PetscFunctionReturn(0);
1454bc6112feSHong Zhang }
1455bc6112feSHong Zhang 
14568928b65cSHong Zhang /* -------------------------------------------------------------------------------------------*/
14578928b65cSHong Zhang #undef __FUNCT__
14588928b65cSHong Zhang #define __FUNCT__ "MatMumpsSetCntl_MUMPS"
14598928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal val)
14608928b65cSHong Zhang {
14618928b65cSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
14628928b65cSHong Zhang 
14638928b65cSHong Zhang   PetscFunctionBegin;
14648928b65cSHong Zhang   mumps->id.CNTL(icntl) = val;
14658928b65cSHong Zhang   PetscFunctionReturn(0);
14668928b65cSHong Zhang }
14678928b65cSHong Zhang 
14688928b65cSHong Zhang #undef __FUNCT__
1469bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetCntl_MUMPS"
1470bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal *val)
1471bc6112feSHong Zhang {
1472bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1473bc6112feSHong Zhang 
1474bc6112feSHong Zhang   PetscFunctionBegin;
1475bc6112feSHong Zhang   *val = mumps->id.CNTL(icntl);
1476bc6112feSHong Zhang   PetscFunctionReturn(0);
1477bc6112feSHong Zhang }
1478bc6112feSHong Zhang 
1479bc6112feSHong Zhang #undef __FUNCT__
14808928b65cSHong Zhang #define __FUNCT__ "MatMumpsSetCntl"
14818928b65cSHong Zhang /*@
14828928b65cSHong Zhang   MatMumpsSetCntl - Set MUMPS parameter CNTL()
14838928b65cSHong Zhang 
14848928b65cSHong Zhang    Logically Collective on Mat
14858928b65cSHong Zhang 
14868928b65cSHong Zhang    Input Parameters:
14878928b65cSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
14888928b65cSHong Zhang .  icntl - index of MUMPS parameter array CNTL()
14898928b65cSHong Zhang -  val - value of MUMPS CNTL(icntl)
14908928b65cSHong Zhang 
14918928b65cSHong Zhang   Options Database:
14928928b65cSHong Zhang .   -mat_mumps_cntl_<icntl> <val>
14938928b65cSHong Zhang 
14948928b65cSHong Zhang    Level: beginner
14958928b65cSHong Zhang 
14968928b65cSHong Zhang    References: MUMPS Users' Guide
14978928b65cSHong Zhang 
14988928b65cSHong Zhang .seealso: MatGetFactor()
14998928b65cSHong Zhang @*/
15008928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl(Mat F,PetscInt icntl,PetscReal val)
15018928b65cSHong Zhang {
15028928b65cSHong Zhang   PetscErrorCode ierr;
15038928b65cSHong Zhang 
15048928b65cSHong Zhang   PetscFunctionBegin;
15058928b65cSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
1506bc6112feSHong Zhang   PetscValidLogicalCollectiveReal(F,val,3);
15078928b65cSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsSetCntl_C",(Mat,PetscInt,PetscReal),(F,icntl,val));CHKERRQ(ierr);
15088928b65cSHong Zhang   PetscFunctionReturn(0);
15098928b65cSHong Zhang }
15108928b65cSHong Zhang 
1511bc6112feSHong Zhang #undef __FUNCT__
1512bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetCntl"
1513a21f80fcSHong Zhang /*@
1514a21f80fcSHong Zhang   MatMumpsGetCntl - Get MUMPS parameter CNTL()
1515a21f80fcSHong Zhang 
1516a21f80fcSHong Zhang    Logically Collective on Mat
1517a21f80fcSHong Zhang 
1518a21f80fcSHong Zhang    Input Parameters:
1519a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1520a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array CNTL()
1521a21f80fcSHong Zhang 
1522a21f80fcSHong Zhang   Output Parameter:
1523a21f80fcSHong Zhang .  val - value of MUMPS CNTL(icntl)
1524a21f80fcSHong Zhang 
1525a21f80fcSHong Zhang    Level: beginner
1526a21f80fcSHong Zhang 
1527a21f80fcSHong Zhang    References: MUMPS Users' Guide
1528a21f80fcSHong Zhang 
1529a21f80fcSHong Zhang .seealso: MatGetFactor()
1530a21f80fcSHong Zhang @*/
1531bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl(Mat F,PetscInt icntl,PetscReal *val)
1532bc6112feSHong Zhang {
1533bc6112feSHong Zhang   PetscErrorCode ierr;
1534bc6112feSHong Zhang 
1535bc6112feSHong Zhang   PetscFunctionBegin;
1536bc6112feSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
1537bc6112feSHong Zhang   PetscValidRealPointer(val,3);
1538bc6112feSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetCntl_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
1539bc6112feSHong Zhang   PetscFunctionReturn(0);
1540bc6112feSHong Zhang }
1541bc6112feSHong Zhang 
1542bc6112feSHong Zhang #undef __FUNCT__
1543ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfo_MUMPS"
1544ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F,PetscInt icntl,PetscInt *info)
1545bc6112feSHong Zhang {
1546bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1547bc6112feSHong Zhang 
1548bc6112feSHong Zhang   PetscFunctionBegin;
1549bc6112feSHong Zhang   *info = mumps->id.INFO(icntl);
1550bc6112feSHong Zhang   PetscFunctionReturn(0);
1551bc6112feSHong Zhang }
1552bc6112feSHong Zhang 
1553bc6112feSHong Zhang #undef __FUNCT__
1554ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfog_MUMPS"
1555ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F,PetscInt icntl,PetscInt *infog)
1556bc6112feSHong Zhang {
1557bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1558bc6112feSHong Zhang 
1559bc6112feSHong Zhang   PetscFunctionBegin;
1560bc6112feSHong Zhang   *infog = mumps->id.INFOG(icntl);
1561bc6112feSHong Zhang   PetscFunctionReturn(0);
1562bc6112feSHong Zhang }
1563bc6112feSHong Zhang 
1564bc6112feSHong Zhang #undef __FUNCT__
1565ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfo_MUMPS"
1566ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfo)
1567bc6112feSHong Zhang {
1568bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1569bc6112feSHong Zhang 
1570bc6112feSHong Zhang   PetscFunctionBegin;
1571bc6112feSHong Zhang   *rinfo = mumps->id.RINFO(icntl);
1572bc6112feSHong Zhang   PetscFunctionReturn(0);
1573bc6112feSHong Zhang }
1574bc6112feSHong Zhang 
1575bc6112feSHong Zhang #undef __FUNCT__
1576ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfog_MUMPS"
1577ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfog)
1578bc6112feSHong Zhang {
1579bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1580bc6112feSHong Zhang 
1581bc6112feSHong Zhang   PetscFunctionBegin;
1582bc6112feSHong Zhang   *rinfog = mumps->id.RINFOG(icntl);
1583bc6112feSHong Zhang   PetscFunctionReturn(0);
1584bc6112feSHong Zhang }
1585bc6112feSHong Zhang 
1586bc6112feSHong Zhang #undef __FUNCT__
1587ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfo"
1588a21f80fcSHong Zhang /*@
1589a21f80fcSHong Zhang   MatMumpsGetInfo - Get MUMPS parameter INFO()
1590a21f80fcSHong Zhang 
1591a21f80fcSHong Zhang    Logically Collective on Mat
1592a21f80fcSHong Zhang 
1593a21f80fcSHong Zhang    Input Parameters:
1594a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1595a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array INFO()
1596a21f80fcSHong Zhang 
1597a21f80fcSHong Zhang   Output Parameter:
1598a21f80fcSHong Zhang .  ival - value of MUMPS INFO(icntl)
1599a21f80fcSHong Zhang 
1600a21f80fcSHong Zhang    Level: beginner
1601a21f80fcSHong Zhang 
1602a21f80fcSHong Zhang    References: MUMPS Users' Guide
1603a21f80fcSHong Zhang 
1604a21f80fcSHong Zhang .seealso: MatGetFactor()
1605a21f80fcSHong Zhang @*/
1606ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo(Mat F,PetscInt icntl,PetscInt *ival)
1607bc6112feSHong Zhang {
1608bc6112feSHong Zhang   PetscErrorCode ierr;
1609bc6112feSHong Zhang 
1610bc6112feSHong Zhang   PetscFunctionBegin;
1611ca810319SHong Zhang   PetscValidIntPointer(ival,3);
1612ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetInfo_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
1613bc6112feSHong Zhang   PetscFunctionReturn(0);
1614bc6112feSHong Zhang }
1615bc6112feSHong Zhang 
1616bc6112feSHong Zhang #undef __FUNCT__
1617ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfog"
1618a21f80fcSHong Zhang /*@
1619a21f80fcSHong Zhang   MatMumpsGetInfog - Get MUMPS parameter INFOG()
1620a21f80fcSHong Zhang 
1621a21f80fcSHong Zhang    Logically Collective on Mat
1622a21f80fcSHong Zhang 
1623a21f80fcSHong Zhang    Input Parameters:
1624a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1625a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array INFOG()
1626a21f80fcSHong Zhang 
1627a21f80fcSHong Zhang   Output Parameter:
1628a21f80fcSHong Zhang .  ival - value of MUMPS INFOG(icntl)
1629a21f80fcSHong Zhang 
1630a21f80fcSHong Zhang    Level: beginner
1631a21f80fcSHong Zhang 
1632a21f80fcSHong Zhang    References: MUMPS Users' Guide
1633a21f80fcSHong Zhang 
1634a21f80fcSHong Zhang .seealso: MatGetFactor()
1635a21f80fcSHong Zhang @*/
1636ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog(Mat F,PetscInt icntl,PetscInt *ival)
1637bc6112feSHong Zhang {
1638bc6112feSHong Zhang   PetscErrorCode ierr;
1639bc6112feSHong Zhang 
1640bc6112feSHong Zhang   PetscFunctionBegin;
1641ca810319SHong Zhang   PetscValidIntPointer(ival,3);
1642ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetInfog_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
1643bc6112feSHong Zhang   PetscFunctionReturn(0);
1644bc6112feSHong Zhang }
1645bc6112feSHong Zhang 
1646bc6112feSHong Zhang #undef __FUNCT__
1647ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfo"
1648a21f80fcSHong Zhang /*@
1649a21f80fcSHong Zhang   MatMumpsGetRinfo - Get MUMPS parameter RINFO()
1650a21f80fcSHong Zhang 
1651a21f80fcSHong Zhang    Logically Collective on Mat
1652a21f80fcSHong Zhang 
1653a21f80fcSHong Zhang    Input Parameters:
1654a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1655a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array RINFO()
1656a21f80fcSHong Zhang 
1657a21f80fcSHong Zhang   Output Parameter:
1658a21f80fcSHong Zhang .  val - value of MUMPS RINFO(icntl)
1659a21f80fcSHong Zhang 
1660a21f80fcSHong Zhang    Level: beginner
1661a21f80fcSHong Zhang 
1662a21f80fcSHong Zhang    References: MUMPS Users' Guide
1663a21f80fcSHong Zhang 
1664a21f80fcSHong Zhang .seealso: MatGetFactor()
1665a21f80fcSHong Zhang @*/
1666ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo(Mat F,PetscInt icntl,PetscReal *val)
1667bc6112feSHong Zhang {
1668bc6112feSHong Zhang   PetscErrorCode ierr;
1669bc6112feSHong Zhang 
1670bc6112feSHong Zhang   PetscFunctionBegin;
1671bc6112feSHong Zhang   PetscValidRealPointer(val,3);
1672ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetRinfo_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
1673bc6112feSHong Zhang   PetscFunctionReturn(0);
1674bc6112feSHong Zhang }
1675bc6112feSHong Zhang 
1676bc6112feSHong Zhang #undef __FUNCT__
1677ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfog"
1678a21f80fcSHong Zhang /*@
1679a21f80fcSHong Zhang   MatMumpsGetRinfog - Get MUMPS parameter RINFOG()
1680a21f80fcSHong Zhang 
1681a21f80fcSHong Zhang    Logically Collective on Mat
1682a21f80fcSHong Zhang 
1683a21f80fcSHong Zhang    Input Parameters:
1684a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1685a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array RINFOG()
1686a21f80fcSHong Zhang 
1687a21f80fcSHong Zhang   Output Parameter:
1688a21f80fcSHong Zhang .  val - value of MUMPS RINFOG(icntl)
1689a21f80fcSHong Zhang 
1690a21f80fcSHong Zhang    Level: beginner
1691a21f80fcSHong Zhang 
1692a21f80fcSHong Zhang    References: MUMPS Users' Guide
1693a21f80fcSHong Zhang 
1694a21f80fcSHong Zhang .seealso: MatGetFactor()
1695a21f80fcSHong Zhang @*/
1696ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog(Mat F,PetscInt icntl,PetscReal *val)
1697bc6112feSHong Zhang {
1698bc6112feSHong Zhang   PetscErrorCode ierr;
1699bc6112feSHong Zhang 
1700bc6112feSHong Zhang   PetscFunctionBegin;
1701bc6112feSHong Zhang   PetscValidRealPointer(val,3);
1702ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetRinfog_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
1703bc6112feSHong Zhang   PetscFunctionReturn(0);
1704bc6112feSHong Zhang }
1705bc6112feSHong Zhang 
170624b6179bSKris Buschelman /*MC
17072692d6eeSBarry Smith   MATSOLVERMUMPS -  A matrix type providing direct solvers (LU and Cholesky) for
170824b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
170924b6179bSKris Buschelman 
171041c8de11SBarry Smith   Works with MATAIJ and MATSBAIJ matrices
171124b6179bSKris Buschelman 
171224b6179bSKris Buschelman   Options Database Keys:
17134e34a73bSHong Zhang +  -mat_mumps_icntl_1 <6>: ICNTL(1): output stream for error messages (None)
17144e34a73bSHong Zhang .  -mat_mumps_icntl_2 <0>: ICNTL(2): output stream for diagnostic printing, statistics, and warning (None)
17154e34a73bSHong Zhang .  -mat_mumps_icntl_3 <0>: ICNTL(3): output stream for global information, collected on the host (None)
17164e34a73bSHong Zhang .  -mat_mumps_icntl_4 <0>: ICNTL(4): level of printing (0 to 4) (None)
17174e34a73bSHong Zhang .  -mat_mumps_icntl_6 <7>: ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7) (None)
17184e34a73bSHong 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)
17194e34a73bSHong Zhang .  -mat_mumps_icntl_8 <77>: ICNTL(8): scaling strategy (-2 to 8 or 77) (None)
17204e34a73bSHong Zhang .  -mat_mumps_icntl_10 <0>: ICNTL(10): max num of refinements (None)
17214e34a73bSHong Zhang .  -mat_mumps_icntl_11 <0>: ICNTL(11): statistics related to an error analysis (via -ksp_view) (None)
17224e34a73bSHong Zhang .  -mat_mumps_icntl_12 <1>: ICNTL(12): an ordering strategy for symmetric matrices (0 to 3) (None)
17234e34a73bSHong Zhang .  -mat_mumps_icntl_13 <0>: ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting (None)
17244e34a73bSHong Zhang .  -mat_mumps_icntl_14 <20>: ICNTL(14): percentage increase in the estimated working space (None)
17254e34a73bSHong Zhang .  -mat_mumps_icntl_19 <0>: ICNTL(19): computes the Schur complement (None)
17264e34a73bSHong Zhang .  -mat_mumps_icntl_22 <0>: ICNTL(22): in-core/out-of-core factorization and solve (0 or 1) (None)
17274e34a73bSHong Zhang .  -mat_mumps_icntl_23 <0>: ICNTL(23): max size of the working memory (MB) that can allocate per processor (None)
17284e34a73bSHong Zhang .  -mat_mumps_icntl_24 <0>: ICNTL(24): detection of null pivot rows (0 or 1) (None)
17294e34a73bSHong Zhang .  -mat_mumps_icntl_25 <0>: ICNTL(25): compute a solution of a deficient matrix and a null space basis (None)
17304e34a73bSHong Zhang .  -mat_mumps_icntl_26 <0>: ICNTL(26): drives the solution phase if a Schur complement matrix (None)
17314e34a73bSHong 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)
17324e34a73bSHong Zhang .  -mat_mumps_icntl_29 <0>: ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis (None)
17334e34a73bSHong Zhang .  -mat_mumps_icntl_30 <0>: ICNTL(30): compute user-specified set of entries in inv(A) (None)
17344e34a73bSHong Zhang .  -mat_mumps_icntl_31 <0>: ICNTL(31): indicates which factors may be discarded during factorization (None)
17354e34a73bSHong Zhang .  -mat_mumps_icntl_33 <0>: ICNTL(33): compute determinant (None)
17364e34a73bSHong Zhang .  -mat_mumps_cntl_1 <0.01>: CNTL(1): relative pivoting threshold (None)
17374e34a73bSHong Zhang .  -mat_mumps_cntl_2 <1.49012e-08>: CNTL(2): stopping criterion of refinement (None)
17384e34a73bSHong Zhang .  -mat_mumps_cntl_3 <0>: CNTL(3): absolute pivoting threshold (None)
17394e34a73bSHong Zhang .  -mat_mumps_cntl_4 <-1>: CNTL(4): value for static pivoting (None)
17404e34a73bSHong Zhang -  -mat_mumps_cntl_5 <0>: CNTL(5): fixation for null pivots (None)
174124b6179bSKris Buschelman 
174224b6179bSKris Buschelman   Level: beginner
174324b6179bSKris Buschelman 
174441c8de11SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage
174541c8de11SBarry Smith 
174624b6179bSKris Buschelman M*/
174724b6179bSKris Buschelman 
174835bd34faSBarry Smith #undef __FUNCT__
174935bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_mumps"
1750f7a08781SBarry Smith static PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type)
175135bd34faSBarry Smith {
175235bd34faSBarry Smith   PetscFunctionBegin;
17532692d6eeSBarry Smith   *type = MATSOLVERMUMPS;
175435bd34faSBarry Smith   PetscFunctionReturn(0);
175535bd34faSBarry Smith }
175635bd34faSBarry Smith 
1757bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI AIJ matrices */
17582877fffaSHong Zhang #undef __FUNCT__
1759bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_aij_mumps"
17608cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
17612877fffaSHong Zhang {
17622877fffaSHong Zhang   Mat            B;
17632877fffaSHong Zhang   PetscErrorCode ierr;
17642877fffaSHong Zhang   Mat_MUMPS      *mumps;
1765ace3abfcSBarry Smith   PetscBool      isSeqAIJ;
17662877fffaSHong Zhang 
17672877fffaSHong Zhang   PetscFunctionBegin;
17682877fffaSHong Zhang   /* Create the factorization matrix */
1769251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr);
1770ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
17712877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
17722877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1773bccb9932SShri Abhyankar   if (isSeqAIJ) {
17740298fd71SBarry Smith     ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr);
1775bccb9932SShri Abhyankar   } else {
17760298fd71SBarry Smith     ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr);
1777bccb9932SShri Abhyankar   }
17782877fffaSHong Zhang 
1779b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
17802205254eSKarl Rupp 
17812877fffaSHong Zhang   B->ops->view        = MatView_MUMPS;
178235bd34faSBarry Smith   B->ops->getinfo     = MatGetInfo_MUMPS;
178320be8e61SHong Zhang   B->ops->getdiagonal = MatGetDiagonal_MUMPS;
17842205254eSKarl Rupp 
1785bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
1786bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
1787bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
1788bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
1789bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
1790bc6112feSHong Zhang 
1791ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
1792ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
1793ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
1794ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
1795450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
1796450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
1797d5f3da31SBarry Smith     B->factortype            = MAT_FACTOR_LU;
1798bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
1799bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
1800746480a1SHong Zhang     mumps->sym = 0;
1801dcd589f8SShri Abhyankar   } else {
180267877ebaSShri Abhyankar     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
1803450b117fSShri Abhyankar     B->factortype                  = MAT_FACTOR_CHOLESKY;
1804bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
1805bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
18066fdc2a6dSBarry Smith     if (A->spd_set && A->spd) mumps->sym = 1;
18076fdc2a6dSBarry Smith     else                      mumps->sym = 2;
1808450b117fSShri Abhyankar   }
18092877fffaSHong Zhang 
18102877fffaSHong Zhang   mumps->isAIJ    = PETSC_TRUE;
1811bf0cc555SLisandro Dalcin   mumps->Destroy  = B->ops->destroy;
18122877fffaSHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
18132877fffaSHong Zhang   B->spptr        = (void*)mumps;
18142205254eSKarl Rupp 
1815f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1816746480a1SHong Zhang 
18172877fffaSHong Zhang   *F = B;
18182877fffaSHong Zhang   PetscFunctionReturn(0);
18192877fffaSHong Zhang }
18202877fffaSHong Zhang 
1821bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI SBAIJ matrices */
18222877fffaSHong Zhang #undef __FUNCT__
1823bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_sbaij_mumps"
18248cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
18252877fffaSHong Zhang {
18262877fffaSHong Zhang   Mat            B;
18272877fffaSHong Zhang   PetscErrorCode ierr;
18282877fffaSHong Zhang   Mat_MUMPS      *mumps;
1829ace3abfcSBarry Smith   PetscBool      isSeqSBAIJ;
18302877fffaSHong Zhang 
18312877fffaSHong Zhang   PetscFunctionBegin;
1832ce94432eSBarry Smith   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
1833ce94432eSBarry 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");
1834251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
18352877fffaSHong Zhang   /* Create the factorization matrix */
1836ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
18372877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
18382877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1839b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
1840bccb9932SShri Abhyankar   if (isSeqSBAIJ) {
18410298fd71SBarry Smith     ierr = MatSeqSBAIJSetPreallocation(B,1,0,NULL);CHKERRQ(ierr);
18422205254eSKarl Rupp 
184316ebf90aSShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
1844dcd589f8SShri Abhyankar   } else {
18450298fd71SBarry Smith     ierr = MatMPISBAIJSetPreallocation(B,1,0,NULL,0,NULL);CHKERRQ(ierr);
18462205254eSKarl Rupp 
1847bccb9932SShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
1848bccb9932SShri Abhyankar   }
1849bccb9932SShri Abhyankar 
185067877ebaSShri Abhyankar   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
1851bccb9932SShri Abhyankar   B->ops->view                   = MatView_MUMPS;
185220be8e61SHong Zhang   B->ops->getdiagonal            = MatGetDiagonal_MUMPS;
18532205254eSKarl Rupp 
1854bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
1855b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
1856b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
1857b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
1858b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
1859bc6112feSHong Zhang 
1860ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
1861ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
1862ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
1863ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
18642205254eSKarl Rupp 
1865f4762488SHong Zhang   B->factortype = MAT_FACTOR_CHOLESKY;
18666fdc2a6dSBarry Smith   if (A->spd_set && A->spd) mumps->sym = 1;
18676fdc2a6dSBarry Smith   else                      mumps->sym = 2;
1868a214ac2aSShri Abhyankar 
1869bccb9932SShri Abhyankar   mumps->isAIJ    = PETSC_FALSE;
1870bf0cc555SLisandro Dalcin   mumps->Destroy  = B->ops->destroy;
1871f3c0ef26SHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
18722877fffaSHong Zhang   B->spptr        = (void*)mumps;
18732205254eSKarl Rupp 
1874f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1875746480a1SHong Zhang 
18762877fffaSHong Zhang   *F = B;
18772877fffaSHong Zhang   PetscFunctionReturn(0);
18782877fffaSHong Zhang }
187997969023SHong Zhang 
1880450b117fSShri Abhyankar #undef __FUNCT__
1881bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_baij_mumps"
18828cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
188367877ebaSShri Abhyankar {
188467877ebaSShri Abhyankar   Mat            B;
188567877ebaSShri Abhyankar   PetscErrorCode ierr;
188667877ebaSShri Abhyankar   Mat_MUMPS      *mumps;
1887ace3abfcSBarry Smith   PetscBool      isSeqBAIJ;
188867877ebaSShri Abhyankar 
188967877ebaSShri Abhyankar   PetscFunctionBegin;
189067877ebaSShri Abhyankar   /* Create the factorization matrix */
1891251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr);
1892ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
189367877ebaSShri Abhyankar   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
189467877ebaSShri Abhyankar   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1895bccb9932SShri Abhyankar   if (isSeqBAIJ) {
18960298fd71SBarry Smith     ierr = MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,NULL);CHKERRQ(ierr);
1897bccb9932SShri Abhyankar   } else {
18980298fd71SBarry Smith     ierr = MatMPIBAIJSetPreallocation(B,A->rmap->bs,0,NULL,0,NULL);CHKERRQ(ierr);
1899bccb9932SShri Abhyankar   }
1900450b117fSShri Abhyankar 
1901b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
1902450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
1903450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
1904450b117fSShri Abhyankar     B->factortype            = MAT_FACTOR_LU;
1905bccb9932SShri Abhyankar     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
1906bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
1907746480a1SHong Zhang     mumps->sym = 0;
1908f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n");
1909bccb9932SShri Abhyankar 
1910450b117fSShri Abhyankar   B->ops->view        = MatView_MUMPS;
191120be8e61SHong Zhang   B->ops->getdiagonal = MatGetDiagonal_MUMPS;
19122205254eSKarl Rupp 
1913bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
1914bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
1915bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
1916bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
1917bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
1918bc6112feSHong Zhang 
1919ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
1920ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
1921ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
1922ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
1923450b117fSShri Abhyankar 
1924450b117fSShri Abhyankar   mumps->isAIJ    = PETSC_TRUE;
1925bf0cc555SLisandro Dalcin   mumps->Destroy  = B->ops->destroy;
1926450b117fSShri Abhyankar   B->ops->destroy = MatDestroy_MUMPS;
1927450b117fSShri Abhyankar   B->spptr        = (void*)mumps;
19282205254eSKarl Rupp 
1929f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1930746480a1SHong Zhang 
1931450b117fSShri Abhyankar   *F = B;
1932450b117fSShri Abhyankar   PetscFunctionReturn(0);
1933450b117fSShri Abhyankar }
193442c9c57cSBarry Smith 
193542c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mumps(Mat,MatFactorType,Mat*);
193642c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_baij_mumps(Mat,MatFactorType,Mat*);
193742c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat,MatFactorType,Mat*);
193842c9c57cSBarry Smith 
193942c9c57cSBarry Smith #undef __FUNCT__
194042c9c57cSBarry Smith #define __FUNCT__ "MatSolverPackageRegister_MUMPS"
194129b38603SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_MUMPS(void)
194242c9c57cSBarry Smith {
194342c9c57cSBarry Smith   PetscErrorCode ierr;
194442c9c57cSBarry Smith 
194542c9c57cSBarry Smith   PetscFunctionBegin;
194642c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ,         MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr);
194742c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ,         MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr);
194842c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ,         MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr);
194942c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ,         MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr);
195042c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPISBAIJ,         MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr);
195142c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ,         MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr);
195242c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ,         MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr);
195342c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ,         MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr);
195442c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ,         MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr);
195542c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQSBAIJ,         MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr);
195642c9c57cSBarry Smith   PetscFunctionReturn(0);
195742c9c57cSBarry Smith }
195842c9c57cSBarry Smith 
1959