xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision 71aed81d4d5f7c964fc3958a492411ca2a96c815)
11c2a3de1SBarry Smith 
2397b6df1SKris Buschelman /*
3c2b5dc30SHong Zhang     Provides an interface to the MUMPS sparse solver
4397b6df1SKris Buschelman */
551d5961aSHong Zhang 
6c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> /*I  "petscmat.h"  I*/
7c6db04a5SJed Brown #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
8397b6df1SKris Buschelman 
9397b6df1SKris Buschelman EXTERN_C_BEGIN
10397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
112907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
122907cef9SHong Zhang #include <cmumps_c.h>
132907cef9SHong Zhang #else
14c6db04a5SJed Brown #include <zmumps_c.h>
152907cef9SHong Zhang #endif
162907cef9SHong Zhang #else
172907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
182907cef9SHong Zhang #include <smumps_c.h>
19397b6df1SKris Buschelman #else
20c6db04a5SJed Brown #include <dmumps_c.h>
21397b6df1SKris Buschelman #endif
222907cef9SHong Zhang #endif
23397b6df1SKris Buschelman EXTERN_C_END
24397b6df1SKris Buschelman #define JOB_INIT -1
253d472b54SHong Zhang #define JOB_FACTSYMBOLIC 1
263d472b54SHong Zhang #define JOB_FACTNUMERIC 2
273d472b54SHong Zhang #define JOB_SOLVE 3
28397b6df1SKris Buschelman #define JOB_END -2
293d472b54SHong Zhang 
302907cef9SHong Zhang /* calls to MUMPS */
312907cef9SHong Zhang #if defined(PETSC_USE_COMPLEX)
322907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
332907cef9SHong Zhang #define PetscMUMPS_c cmumps_c
342907cef9SHong Zhang #else
352907cef9SHong Zhang #define PetscMUMPS_c zmumps_c
362907cef9SHong Zhang #endif
372907cef9SHong Zhang #else
382907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
392907cef9SHong Zhang #define PetscMUMPS_c smumps_c
402907cef9SHong Zhang #else
412907cef9SHong Zhang #define PetscMUMPS_c dmumps_c
422907cef9SHong Zhang #endif
432907cef9SHong Zhang #endif
442907cef9SHong Zhang 
453d472b54SHong Zhang 
46397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */
47397b6df1SKris Buschelman #define ICNTL(I) icntl[(I)-1]
48397b6df1SKris Buschelman #define CNTL(I) cntl[(I)-1]
49397b6df1SKris Buschelman #define INFOG(I) infog[(I)-1]
50a7aca84bSHong Zhang #define INFO(I) info[(I)-1]
51397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1]
52adc1d99fSHong Zhang #define RINFO(I) rinfo[(I)-1]
53397b6df1SKris Buschelman 
54397b6df1SKris Buschelman typedef struct {
55397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
562907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
572907cef9SHong Zhang   CMUMPS_STRUC_C id;
582907cef9SHong Zhang #else
59397b6df1SKris Buschelman   ZMUMPS_STRUC_C id;
602907cef9SHong Zhang #endif
612907cef9SHong Zhang #else
622907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
632907cef9SHong Zhang   SMUMPS_STRUC_C id;
64397b6df1SKris Buschelman #else
65397b6df1SKris Buschelman   DMUMPS_STRUC_C id;
66397b6df1SKris Buschelman #endif
672907cef9SHong Zhang #endif
682907cef9SHong Zhang 
69397b6df1SKris Buschelman   MatStructure matstruc;
70c1490034SHong Zhang   PetscMPIInt  myid,size;
71a5e57a09SHong Zhang   PetscInt     *irn,*jcn,nz,sym;
72397b6df1SKris Buschelman   PetscScalar  *val;
73397b6df1SKris Buschelman   MPI_Comm     comm_mumps;
7464e6c443SBarry Smith   PetscBool    isAIJ,CleanUpMUMPS;
75a5e57a09SHong Zhang   PetscInt     ICNTL9_pre;           /* check if ICNTL(9) is changed from previous MatSolve */
76801fbe65SHong Zhang   VecScatter   scat_rhs, scat_sol;   /* used by MatSolve() */
77801fbe65SHong Zhang   Vec          b_seq,x_seq;
782205254eSKarl Rupp 
79bf0cc555SLisandro Dalcin   PetscErrorCode (*Destroy)(Mat);
80bccb9932SShri Abhyankar   PetscErrorCode (*ConvertToTriples)(Mat, int, MatReuse, int*, int**, int**, PetscScalar**);
81f0c56d0fSKris Buschelman } Mat_MUMPS;
82f0c56d0fSKris Buschelman 
8309573ac7SBarry Smith extern PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);
84b24902e0SBarry Smith 
85397b6df1SKris Buschelman /*
86d341cd04SHong Zhang   MatConvertToTriples_A_B - convert Petsc matrix to triples: row[nz], col[nz], val[nz]
87d341cd04SHong Zhang 
88397b6df1SKris Buschelman   input:
8967877ebaSShri Abhyankar     A       - matrix in aij,baij or sbaij (bs=1) format
90397b6df1SKris Buschelman     shift   - 0: C style output triple; 1: Fortran style output triple.
91bccb9932SShri Abhyankar     reuse   - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple
92bccb9932SShri Abhyankar               MAT_REUSE_MATRIX:   only the values in v array are updated
93397b6df1SKris Buschelman   output:
94397b6df1SKris Buschelman     nnz     - dim of r, c, and v (number of local nonzero entries of A)
95397b6df1SKris Buschelman     r, c, v - row and col index, matrix values (matrix triples)
96eb9baa12SBarry Smith 
97eb9baa12SBarry Smith   The returned values r, c, and sometimes v are obtained in a single PetscMalloc(). Then in MatDestroy_MUMPS() it is
98eb9baa12SBarry Smith   freed with PetscFree((mumps->irn);  This is not ideal code, the fact that v is ONLY sometimes part of mumps->irn means
99eb9baa12SBarry Smith   that the PetscMalloc() cannot easily be replaced with a PetscMalloc3().
100eb9baa12SBarry Smith 
101397b6df1SKris Buschelman  */
10216ebf90aSShri Abhyankar 
10316ebf90aSShri Abhyankar #undef __FUNCT__
10416ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqaij"
105bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
106b24902e0SBarry Smith {
107185f6596SHong Zhang   const PetscInt *ai,*aj,*ajj,M=A->rmap->n;
10867877ebaSShri Abhyankar   PetscInt       nz,rnz,i,j;
109dfbe8321SBarry Smith   PetscErrorCode ierr;
110c1490034SHong Zhang   PetscInt       *row,*col;
11116ebf90aSShri Abhyankar   Mat_SeqAIJ     *aa=(Mat_SeqAIJ*)A->data;
112397b6df1SKris Buschelman 
113397b6df1SKris Buschelman   PetscFunctionBegin;
11416ebf90aSShri Abhyankar   *v=aa->a;
115bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
1162205254eSKarl Rupp     nz   = aa->nz;
1172205254eSKarl Rupp     ai   = aa->i;
1182205254eSKarl Rupp     aj   = aa->j;
11916ebf90aSShri Abhyankar     *nnz = nz;
120785e854fSJed Brown     ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr);
121185f6596SHong Zhang     col  = row + nz;
122185f6596SHong Zhang 
12316ebf90aSShri Abhyankar     nz = 0;
12416ebf90aSShri Abhyankar     for (i=0; i<M; i++) {
12516ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
12667877ebaSShri Abhyankar       ajj = aj + ai[i];
12767877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
12867877ebaSShri Abhyankar         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
12916ebf90aSShri Abhyankar       }
13016ebf90aSShri Abhyankar     }
13116ebf90aSShri Abhyankar     *r = row; *c = col;
13216ebf90aSShri Abhyankar   }
13316ebf90aSShri Abhyankar   PetscFunctionReturn(0);
13416ebf90aSShri Abhyankar }
135397b6df1SKris Buschelman 
13616ebf90aSShri Abhyankar #undef __FUNCT__
13767877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqbaij_seqaij"
138bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
13967877ebaSShri Abhyankar {
14067877ebaSShri Abhyankar   Mat_SeqBAIJ    *aa=(Mat_SeqBAIJ*)A->data;
14133d57670SJed Brown   const PetscInt *ai,*aj,*ajj,bs2 = aa->bs2;
14233d57670SJed Brown   PetscInt       bs,M,nz,idx=0,rnz,i,j,k,m;
14367877ebaSShri Abhyankar   PetscErrorCode ierr;
14467877ebaSShri Abhyankar   PetscInt       *row,*col;
14567877ebaSShri Abhyankar 
14667877ebaSShri Abhyankar   PetscFunctionBegin;
14733d57670SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
14833d57670SJed Brown   M = A->rmap->N/bs;
149cf3759fdSShri Abhyankar   *v = aa->a;
150bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
151cf3759fdSShri Abhyankar     ai   = aa->i; aj = aa->j;
15267877ebaSShri Abhyankar     nz   = bs2*aa->nz;
15367877ebaSShri Abhyankar     *nnz = nz;
154785e854fSJed Brown     ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr);
155185f6596SHong Zhang     col  = row + nz;
156185f6596SHong Zhang 
15767877ebaSShri Abhyankar     for (i=0; i<M; i++) {
15867877ebaSShri Abhyankar       ajj = aj + ai[i];
15967877ebaSShri Abhyankar       rnz = ai[i+1] - ai[i];
16067877ebaSShri Abhyankar       for (k=0; k<rnz; k++) {
16167877ebaSShri Abhyankar         for (j=0; j<bs; j++) {
16267877ebaSShri Abhyankar           for (m=0; m<bs; m++) {
16367877ebaSShri Abhyankar             row[idx]   = i*bs + m + shift;
164cf3759fdSShri Abhyankar             col[idx++] = bs*(ajj[k]) + j + shift;
16567877ebaSShri Abhyankar           }
16667877ebaSShri Abhyankar         }
16767877ebaSShri Abhyankar       }
16867877ebaSShri Abhyankar     }
169cf3759fdSShri Abhyankar     *r = row; *c = col;
17067877ebaSShri Abhyankar   }
17167877ebaSShri Abhyankar   PetscFunctionReturn(0);
17267877ebaSShri Abhyankar }
17367877ebaSShri Abhyankar 
17467877ebaSShri Abhyankar #undef __FUNCT__
17516ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqsbaij_seqsbaij"
176bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
17716ebf90aSShri Abhyankar {
17867877ebaSShri Abhyankar   const PetscInt *ai, *aj,*ajj,M=A->rmap->n;
17967877ebaSShri Abhyankar   PetscInt       nz,rnz,i,j;
18016ebf90aSShri Abhyankar   PetscErrorCode ierr;
18116ebf90aSShri Abhyankar   PetscInt       *row,*col;
18216ebf90aSShri Abhyankar   Mat_SeqSBAIJ   *aa=(Mat_SeqSBAIJ*)A->data;
18316ebf90aSShri Abhyankar 
18416ebf90aSShri Abhyankar   PetscFunctionBegin;
185882afa5aSHong Zhang   *v = aa->a;
186bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
1872205254eSKarl Rupp     nz   = aa->nz;
1882205254eSKarl Rupp     ai   = aa->i;
1892205254eSKarl Rupp     aj   = aa->j;
1902205254eSKarl Rupp     *v   = aa->a;
19116ebf90aSShri Abhyankar     *nnz = nz;
192785e854fSJed Brown     ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr);
193185f6596SHong Zhang     col  = row + nz;
194185f6596SHong Zhang 
19516ebf90aSShri Abhyankar     nz = 0;
19616ebf90aSShri Abhyankar     for (i=0; i<M; i++) {
19716ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
19867877ebaSShri Abhyankar       ajj = aj + ai[i];
19967877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
20067877ebaSShri Abhyankar         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
20116ebf90aSShri Abhyankar       }
20216ebf90aSShri Abhyankar     }
20316ebf90aSShri Abhyankar     *r = row; *c = col;
20416ebf90aSShri Abhyankar   }
20516ebf90aSShri Abhyankar   PetscFunctionReturn(0);
20616ebf90aSShri Abhyankar }
20716ebf90aSShri Abhyankar 
20816ebf90aSShri Abhyankar #undef __FUNCT__
20916ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqsbaij"
210bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
21116ebf90aSShri Abhyankar {
21267877ebaSShri Abhyankar   const PetscInt    *ai,*aj,*ajj,*adiag,M=A->rmap->n;
21367877ebaSShri Abhyankar   PetscInt          nz,rnz,i,j;
21467877ebaSShri Abhyankar   const PetscScalar *av,*v1;
21516ebf90aSShri Abhyankar   PetscScalar       *val;
21616ebf90aSShri Abhyankar   PetscErrorCode    ierr;
21716ebf90aSShri Abhyankar   PetscInt          *row,*col;
218829b1710SHong Zhang   Mat_SeqAIJ        *aa=(Mat_SeqAIJ*)A->data;
21916ebf90aSShri Abhyankar 
22016ebf90aSShri Abhyankar   PetscFunctionBegin;
22116ebf90aSShri Abhyankar   ai   =aa->i; aj=aa->j;av=aa->a;
22216ebf90aSShri Abhyankar   adiag=aa->diag;
223bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
224829b1710SHong Zhang     /* count nz in the uppper triangular part of A */
225829b1710SHong Zhang     nz = 0;
226829b1710SHong Zhang     for (i=0; i<M; i++) nz += ai[i+1] - adiag[i];
22716ebf90aSShri Abhyankar     *nnz = nz;
228829b1710SHong Zhang 
229185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
230185f6596SHong Zhang     col  = row + nz;
231185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
232185f6596SHong Zhang 
23316ebf90aSShri Abhyankar     nz = 0;
23416ebf90aSShri Abhyankar     for (i=0; i<M; i++) {
23516ebf90aSShri Abhyankar       rnz = ai[i+1] - adiag[i];
23667877ebaSShri Abhyankar       ajj = aj + adiag[i];
237cf3759fdSShri Abhyankar       v1  = av + adiag[i];
23867877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
23967877ebaSShri Abhyankar         row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
24016ebf90aSShri Abhyankar       }
24116ebf90aSShri Abhyankar     }
24216ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
243397b6df1SKris Buschelman   } else {
24416ebf90aSShri Abhyankar     nz = 0; val = *v;
24516ebf90aSShri Abhyankar     for (i=0; i <M; i++) {
24616ebf90aSShri Abhyankar       rnz = ai[i+1] - adiag[i];
24767877ebaSShri Abhyankar       ajj = aj + adiag[i];
24867877ebaSShri Abhyankar       v1  = av + adiag[i];
24967877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
25067877ebaSShri Abhyankar         val[nz++] = v1[j];
25116ebf90aSShri Abhyankar       }
25216ebf90aSShri Abhyankar     }
25316ebf90aSShri Abhyankar   }
25416ebf90aSShri Abhyankar   PetscFunctionReturn(0);
25516ebf90aSShri Abhyankar }
25616ebf90aSShri Abhyankar 
25716ebf90aSShri Abhyankar #undef __FUNCT__
25816ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpisbaij_mpisbaij"
259bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
26016ebf90aSShri Abhyankar {
26116ebf90aSShri Abhyankar   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
26216ebf90aSShri Abhyankar   PetscErrorCode    ierr;
26316ebf90aSShri Abhyankar   PetscInt          rstart,nz,i,j,jj,irow,countA,countB;
26416ebf90aSShri Abhyankar   PetscInt          *row,*col;
26516ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
26616ebf90aSShri Abhyankar   PetscScalar       *val;
267397b6df1SKris Buschelman   Mat_MPISBAIJ      *mat = (Mat_MPISBAIJ*)A->data;
268397b6df1SKris Buschelman   Mat_SeqSBAIJ      *aa  = (Mat_SeqSBAIJ*)(mat->A)->data;
269397b6df1SKris Buschelman   Mat_SeqBAIJ       *bb  = (Mat_SeqBAIJ*)(mat->B)->data;
27016ebf90aSShri Abhyankar 
27116ebf90aSShri Abhyankar   PetscFunctionBegin;
272d0f46423SBarry Smith   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
273397b6df1SKris Buschelman   av=aa->a; bv=bb->a;
274397b6df1SKris Buschelman 
2752205254eSKarl Rupp   garray = mat->garray;
2762205254eSKarl Rupp 
277bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
27816ebf90aSShri Abhyankar     nz   = aa->nz + bb->nz;
27916ebf90aSShri Abhyankar     *nnz = nz;
280185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
281185f6596SHong Zhang     col  = row + nz;
282185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
283185f6596SHong Zhang 
284397b6df1SKris Buschelman     *r = row; *c = col; *v = val;
285397b6df1SKris Buschelman   } else {
286397b6df1SKris Buschelman     row = *r; col = *c; val = *v;
287397b6df1SKris Buschelman   }
288397b6df1SKris Buschelman 
289028e57e8SHong Zhang   jj = 0; irow = rstart;
290397b6df1SKris Buschelman   for (i=0; i<m; i++) {
291397b6df1SKris Buschelman     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
292397b6df1SKris Buschelman     countA = ai[i+1] - ai[i];
293397b6df1SKris Buschelman     countB = bi[i+1] - bi[i];
294397b6df1SKris Buschelman     bjj    = bj + bi[i];
29516ebf90aSShri Abhyankar     v1     = av + ai[i];
29616ebf90aSShri Abhyankar     v2     = bv + bi[i];
297397b6df1SKris Buschelman 
298397b6df1SKris Buschelman     /* A-part */
299397b6df1SKris Buschelman     for (j=0; j<countA; j++) {
300bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
301397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
302397b6df1SKris Buschelman       }
30316ebf90aSShri Abhyankar       val[jj++] = v1[j];
304397b6df1SKris Buschelman     }
30516ebf90aSShri Abhyankar 
30616ebf90aSShri Abhyankar     /* B-part */
30716ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
308bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
309397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
310397b6df1SKris Buschelman       }
31116ebf90aSShri Abhyankar       val[jj++] = v2[j];
31216ebf90aSShri Abhyankar     }
31316ebf90aSShri Abhyankar     irow++;
31416ebf90aSShri Abhyankar   }
31516ebf90aSShri Abhyankar   PetscFunctionReturn(0);
31616ebf90aSShri Abhyankar }
31716ebf90aSShri Abhyankar 
31816ebf90aSShri Abhyankar #undef __FUNCT__
31916ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpiaij"
320bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
32116ebf90aSShri Abhyankar {
32216ebf90aSShri Abhyankar   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
32316ebf90aSShri Abhyankar   PetscErrorCode    ierr;
32416ebf90aSShri Abhyankar   PetscInt          rstart,nz,i,j,jj,irow,countA,countB;
32516ebf90aSShri Abhyankar   PetscInt          *row,*col;
32616ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
32716ebf90aSShri Abhyankar   PetscScalar       *val;
32816ebf90aSShri Abhyankar   Mat_MPIAIJ        *mat = (Mat_MPIAIJ*)A->data;
32916ebf90aSShri Abhyankar   Mat_SeqAIJ        *aa  = (Mat_SeqAIJ*)(mat->A)->data;
33016ebf90aSShri Abhyankar   Mat_SeqAIJ        *bb  = (Mat_SeqAIJ*)(mat->B)->data;
33116ebf90aSShri Abhyankar 
33216ebf90aSShri Abhyankar   PetscFunctionBegin;
33316ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
33416ebf90aSShri Abhyankar   av=aa->a; bv=bb->a;
33516ebf90aSShri Abhyankar 
3362205254eSKarl Rupp   garray = mat->garray;
3372205254eSKarl Rupp 
338bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
33916ebf90aSShri Abhyankar     nz   = aa->nz + bb->nz;
34016ebf90aSShri Abhyankar     *nnz = nz;
341185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
342185f6596SHong Zhang     col  = row + nz;
343185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
344185f6596SHong Zhang 
34516ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
34616ebf90aSShri Abhyankar   } else {
34716ebf90aSShri Abhyankar     row = *r; col = *c; val = *v;
34816ebf90aSShri Abhyankar   }
34916ebf90aSShri Abhyankar 
35016ebf90aSShri Abhyankar   jj = 0; irow = rstart;
35116ebf90aSShri Abhyankar   for (i=0; i<m; i++) {
35216ebf90aSShri Abhyankar     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
35316ebf90aSShri Abhyankar     countA = ai[i+1] - ai[i];
35416ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
35516ebf90aSShri Abhyankar     bjj    = bj + bi[i];
35616ebf90aSShri Abhyankar     v1     = av + ai[i];
35716ebf90aSShri Abhyankar     v2     = bv + bi[i];
35816ebf90aSShri Abhyankar 
35916ebf90aSShri Abhyankar     /* A-part */
36016ebf90aSShri Abhyankar     for (j=0; j<countA; j++) {
361bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
36216ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
36316ebf90aSShri Abhyankar       }
36416ebf90aSShri Abhyankar       val[jj++] = v1[j];
36516ebf90aSShri Abhyankar     }
36616ebf90aSShri Abhyankar 
36716ebf90aSShri Abhyankar     /* B-part */
36816ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
369bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
37016ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
37116ebf90aSShri Abhyankar       }
37216ebf90aSShri Abhyankar       val[jj++] = v2[j];
37316ebf90aSShri Abhyankar     }
37416ebf90aSShri Abhyankar     irow++;
37516ebf90aSShri Abhyankar   }
37616ebf90aSShri Abhyankar   PetscFunctionReturn(0);
37716ebf90aSShri Abhyankar }
37816ebf90aSShri Abhyankar 
37916ebf90aSShri Abhyankar #undef __FUNCT__
38067877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpibaij_mpiaij"
381bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
38267877ebaSShri Abhyankar {
38367877ebaSShri Abhyankar   Mat_MPIBAIJ       *mat    = (Mat_MPIBAIJ*)A->data;
38467877ebaSShri Abhyankar   Mat_SeqBAIJ       *aa     = (Mat_SeqBAIJ*)(mat->A)->data;
38567877ebaSShri Abhyankar   Mat_SeqBAIJ       *bb     = (Mat_SeqBAIJ*)(mat->B)->data;
38667877ebaSShri Abhyankar   const PetscInt    *ai     = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj;
387d985c460SShri Abhyankar   const PetscInt    *garray = mat->garray,mbs=mat->mbs,rstart=A->rmap->rstart;
38833d57670SJed Brown   const PetscInt    bs2=mat->bs2;
38967877ebaSShri Abhyankar   PetscErrorCode    ierr;
39033d57670SJed Brown   PetscInt          bs,nz,i,j,k,n,jj,irow,countA,countB,idx;
39167877ebaSShri Abhyankar   PetscInt          *row,*col;
39267877ebaSShri Abhyankar   const PetscScalar *av=aa->a, *bv=bb->a,*v1,*v2;
39367877ebaSShri Abhyankar   PetscScalar       *val;
39467877ebaSShri Abhyankar 
39567877ebaSShri Abhyankar   PetscFunctionBegin;
39633d57670SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
397bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
39867877ebaSShri Abhyankar     nz   = bs2*(aa->nz + bb->nz);
39967877ebaSShri Abhyankar     *nnz = nz;
400185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
401185f6596SHong Zhang     col  = row + nz;
402185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
403185f6596SHong Zhang 
40467877ebaSShri Abhyankar     *r = row; *c = col; *v = val;
40567877ebaSShri Abhyankar   } else {
40667877ebaSShri Abhyankar     row = *r; col = *c; val = *v;
40767877ebaSShri Abhyankar   }
40867877ebaSShri Abhyankar 
409d985c460SShri Abhyankar   jj = 0; irow = rstart;
41067877ebaSShri Abhyankar   for (i=0; i<mbs; i++) {
41167877ebaSShri Abhyankar     countA = ai[i+1] - ai[i];
41267877ebaSShri Abhyankar     countB = bi[i+1] - bi[i];
41367877ebaSShri Abhyankar     ajj    = aj + ai[i];
41467877ebaSShri Abhyankar     bjj    = bj + bi[i];
41567877ebaSShri Abhyankar     v1     = av + bs2*ai[i];
41667877ebaSShri Abhyankar     v2     = bv + bs2*bi[i];
41767877ebaSShri Abhyankar 
41867877ebaSShri Abhyankar     idx = 0;
41967877ebaSShri Abhyankar     /* A-part */
42067877ebaSShri Abhyankar     for (k=0; k<countA; k++) {
42167877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
42267877ebaSShri Abhyankar         for (n=0; n<bs; n++) {
423bccb9932SShri Abhyankar           if (reuse == MAT_INITIAL_MATRIX) {
424d985c460SShri Abhyankar             row[jj] = irow + n + shift;
425d985c460SShri Abhyankar             col[jj] = rstart + bs*ajj[k] + j + shift;
42667877ebaSShri Abhyankar           }
42767877ebaSShri Abhyankar           val[jj++] = v1[idx++];
42867877ebaSShri Abhyankar         }
42967877ebaSShri Abhyankar       }
43067877ebaSShri Abhyankar     }
43167877ebaSShri Abhyankar 
43267877ebaSShri Abhyankar     idx = 0;
43367877ebaSShri Abhyankar     /* B-part */
43467877ebaSShri Abhyankar     for (k=0; k<countB; k++) {
43567877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
43667877ebaSShri Abhyankar         for (n=0; n<bs; n++) {
437bccb9932SShri Abhyankar           if (reuse == MAT_INITIAL_MATRIX) {
438d985c460SShri Abhyankar             row[jj] = irow + n + shift;
439d985c460SShri Abhyankar             col[jj] = bs*garray[bjj[k]] + j + shift;
44067877ebaSShri Abhyankar           }
441d985c460SShri Abhyankar           val[jj++] = v2[idx++];
44267877ebaSShri Abhyankar         }
44367877ebaSShri Abhyankar       }
44467877ebaSShri Abhyankar     }
445d985c460SShri Abhyankar     irow += bs;
44667877ebaSShri Abhyankar   }
44767877ebaSShri Abhyankar   PetscFunctionReturn(0);
44867877ebaSShri Abhyankar }
44967877ebaSShri Abhyankar 
45067877ebaSShri Abhyankar #undef __FUNCT__
45116ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpisbaij"
452bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
45316ebf90aSShri Abhyankar {
45416ebf90aSShri Abhyankar   const PetscInt    *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
45516ebf90aSShri Abhyankar   PetscErrorCode    ierr;
456e0bace9bSHong Zhang   PetscInt          rstart,nz,nza,nzb,i,j,jj,irow,countA,countB;
45716ebf90aSShri Abhyankar   PetscInt          *row,*col;
45816ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
45916ebf90aSShri Abhyankar   PetscScalar       *val;
46016ebf90aSShri Abhyankar   Mat_MPIAIJ        *mat =  (Mat_MPIAIJ*)A->data;
46116ebf90aSShri Abhyankar   Mat_SeqAIJ        *aa  =(Mat_SeqAIJ*)(mat->A)->data;
46216ebf90aSShri Abhyankar   Mat_SeqAIJ        *bb  =(Mat_SeqAIJ*)(mat->B)->data;
46316ebf90aSShri Abhyankar 
46416ebf90aSShri Abhyankar   PetscFunctionBegin;
46516ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j; adiag=aa->diag;
46616ebf90aSShri Abhyankar   bi=bb->i; bj=bb->j; garray = mat->garray;
46716ebf90aSShri Abhyankar   av=aa->a; bv=bb->a;
4682205254eSKarl Rupp 
46916ebf90aSShri Abhyankar   rstart = A->rmap->rstart;
47016ebf90aSShri Abhyankar 
471bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
472e0bace9bSHong Zhang     nza = 0;    /* num of upper triangular entries in mat->A, including diagonals */
473e0bace9bSHong Zhang     nzb = 0;    /* num of upper triangular entries in mat->B */
47416ebf90aSShri Abhyankar     for (i=0; i<m; i++) {
475e0bace9bSHong Zhang       nza   += (ai[i+1] - adiag[i]);
47616ebf90aSShri Abhyankar       countB = bi[i+1] - bi[i];
47716ebf90aSShri Abhyankar       bjj    = bj + bi[i];
478e0bace9bSHong Zhang       for (j=0; j<countB; j++) {
479e0bace9bSHong Zhang         if (garray[bjj[j]] > rstart) nzb++;
480e0bace9bSHong Zhang       }
481e0bace9bSHong Zhang     }
48216ebf90aSShri Abhyankar 
483e0bace9bSHong Zhang     nz   = nza + nzb; /* total nz of upper triangular part of mat */
48416ebf90aSShri Abhyankar     *nnz = nz;
485185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
486185f6596SHong Zhang     col  = row + nz;
487185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
488185f6596SHong Zhang 
48916ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
49016ebf90aSShri Abhyankar   } else {
49116ebf90aSShri Abhyankar     row = *r; col = *c; val = *v;
49216ebf90aSShri Abhyankar   }
49316ebf90aSShri Abhyankar 
49416ebf90aSShri Abhyankar   jj = 0; irow = rstart;
49516ebf90aSShri Abhyankar   for (i=0; i<m; i++) {
49616ebf90aSShri Abhyankar     ajj    = aj + adiag[i];                 /* ptr to the beginning of the diagonal of this row */
49716ebf90aSShri Abhyankar     v1     = av + adiag[i];
49816ebf90aSShri Abhyankar     countA = ai[i+1] - adiag[i];
49916ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
50016ebf90aSShri Abhyankar     bjj    = bj + bi[i];
50116ebf90aSShri Abhyankar     v2     = bv + bi[i];
50216ebf90aSShri Abhyankar 
50316ebf90aSShri Abhyankar     /* A-part */
50416ebf90aSShri Abhyankar     for (j=0; j<countA; j++) {
505bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
50616ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
50716ebf90aSShri Abhyankar       }
50816ebf90aSShri Abhyankar       val[jj++] = v1[j];
50916ebf90aSShri Abhyankar     }
51016ebf90aSShri Abhyankar 
51116ebf90aSShri Abhyankar     /* B-part */
51216ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
51316ebf90aSShri Abhyankar       if (garray[bjj[j]] > rstart) {
514bccb9932SShri Abhyankar         if (reuse == MAT_INITIAL_MATRIX) {
51516ebf90aSShri Abhyankar           row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
51616ebf90aSShri Abhyankar         }
51716ebf90aSShri Abhyankar         val[jj++] = v2[j];
51816ebf90aSShri Abhyankar       }
519397b6df1SKris Buschelman     }
520397b6df1SKris Buschelman     irow++;
521397b6df1SKris Buschelman   }
522397b6df1SKris Buschelman   PetscFunctionReturn(0);
523397b6df1SKris Buschelman }
524397b6df1SKris Buschelman 
525397b6df1SKris Buschelman #undef __FUNCT__
52620be8e61SHong Zhang #define __FUNCT__ "MatGetDiagonal_MUMPS"
52720be8e61SHong Zhang PetscErrorCode MatGetDiagonal_MUMPS(Mat A,Vec v)
52820be8e61SHong Zhang {
52920be8e61SHong Zhang   PetscFunctionBegin;
53020be8e61SHong Zhang   SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type: MUMPS factor");
53120be8e61SHong Zhang   PetscFunctionReturn(0);
53220be8e61SHong Zhang }
53320be8e61SHong Zhang 
53420be8e61SHong Zhang #undef __FUNCT__
5353924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS"
536dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A)
537dfbe8321SBarry Smith {
538a5e57a09SHong Zhang   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
539dfbe8321SBarry Smith   PetscErrorCode ierr;
540b24902e0SBarry Smith 
541397b6df1SKris Buschelman   PetscFunctionBegin;
542a5e57a09SHong Zhang   if (mumps->CleanUpMUMPS) {
543397b6df1SKris Buschelman     /* Terminate instance, deallocate memories */
544a5e57a09SHong Zhang     ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr);
545a5e57a09SHong Zhang     ierr = VecScatterDestroy(&mumps->scat_rhs);CHKERRQ(ierr);
546a5e57a09SHong Zhang     ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
547801fbe65SHong Zhang     ierr = VecDestroy(&mumps->b_seq);CHKERRQ(ierr);
548a5e57a09SHong Zhang     ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr);
549a5e57a09SHong Zhang     ierr = PetscFree(mumps->id.perm_in);CHKERRQ(ierr);
550a5e57a09SHong Zhang     ierr = PetscFree(mumps->irn);CHKERRQ(ierr);
5512205254eSKarl Rupp 
552a5e57a09SHong Zhang     mumps->id.job = JOB_END;
553a5e57a09SHong Zhang     PetscMUMPS_c(&mumps->id);
554a5e57a09SHong Zhang     ierr = MPI_Comm_free(&(mumps->comm_mumps));CHKERRQ(ierr);
555397b6df1SKris Buschelman   }
556a5e57a09SHong Zhang   if (mumps->Destroy) {
557a5e57a09SHong Zhang     ierr = (mumps->Destroy)(A);CHKERRQ(ierr);
558bf0cc555SLisandro Dalcin   }
559bf0cc555SLisandro Dalcin   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
560bf0cc555SLisandro Dalcin 
56197969023SHong Zhang   /* clear composed functions */
562bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr);
563bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetIcntl_C",NULL);CHKERRQ(ierr);
564bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetIcntl_C",NULL);CHKERRQ(ierr);
565bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetCntl_C",NULL);CHKERRQ(ierr);
566bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetCntl_C",NULL);CHKERRQ(ierr);
567bc6112feSHong Zhang 
568ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfo_C",NULL);CHKERRQ(ierr);
569ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfog_C",NULL);CHKERRQ(ierr);
570ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfo_C",NULL);CHKERRQ(ierr);
571ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfog_C",NULL);CHKERRQ(ierr);
572397b6df1SKris Buschelman   PetscFunctionReturn(0);
573397b6df1SKris Buschelman }
574397b6df1SKris Buschelman 
575397b6df1SKris Buschelman #undef __FUNCT__
576f6c57405SHong Zhang #define __FUNCT__ "MatSolve_MUMPS"
577b24902e0SBarry Smith PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x)
578b24902e0SBarry Smith {
579a5e57a09SHong Zhang   Mat_MUMPS        *mumps=(Mat_MUMPS*)A->spptr;
580d54de34fSKris Buschelman   PetscScalar      *array;
58167877ebaSShri Abhyankar   Vec              b_seq;
582329ec9b3SHong Zhang   IS               is_iden,is_petsc;
583dfbe8321SBarry Smith   PetscErrorCode   ierr;
584329ec9b3SHong Zhang   PetscInt         i;
585883f2eb9SBarry Smith   static PetscBool cite1 = PETSC_FALSE,cite2 = PETSC_FALSE;
586397b6df1SKris Buschelman 
587397b6df1SKris Buschelman   PetscFunctionBegin;
588883f2eb9SBarry 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);
589883f2eb9SBarry 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);
590a5e57a09SHong Zhang   mumps->id.nrhs = 1;
591a5e57a09SHong Zhang   b_seq          = mumps->b_seq;
592a5e57a09SHong Zhang   if (mumps->size > 1) {
593329ec9b3SHong Zhang     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
594a5e57a09SHong Zhang     ierr = VecScatterBegin(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
595a5e57a09SHong Zhang     ierr = VecScatterEnd(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
596a5e57a09SHong Zhang     if (!mumps->myid) {ierr = VecGetArray(b_seq,&array);CHKERRQ(ierr);}
597397b6df1SKris Buschelman   } else {  /* size == 1 */
598397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
599397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
600397b6df1SKris Buschelman   }
601a5e57a09SHong Zhang   if (!mumps->myid) { /* define rhs on the host */
602a5e57a09SHong Zhang     mumps->id.nrhs = 1;
603397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
6042907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
605a5e57a09SHong Zhang     mumps->id.rhs = (mumps_complex*)array;
6062907cef9SHong Zhang #else
607a5e57a09SHong Zhang     mumps->id.rhs = (mumps_double_complex*)array;
6082907cef9SHong Zhang #endif
609397b6df1SKris Buschelman #else
610a5e57a09SHong Zhang     mumps->id.rhs = array;
611397b6df1SKris Buschelman #endif
612397b6df1SKris Buschelman   }
613397b6df1SKris Buschelman 
614397b6df1SKris Buschelman   /* solve phase */
615329ec9b3SHong Zhang   /*-------------*/
616a5e57a09SHong Zhang   mumps->id.job = JOB_SOLVE;
617a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
618a5e57a09SHong 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));
619397b6df1SKris Buschelman 
620a5e57a09SHong Zhang   if (mumps->size > 1) { /* convert mumps distributed solution to petsc mpi x */
621a5e57a09SHong Zhang     if (mumps->scat_sol && mumps->ICNTL9_pre != mumps->id.ICNTL(9)) {
622a5e57a09SHong Zhang       /* when id.ICNTL(9) changes, the contents of lsol_loc may change (not its size, lsol_loc), recreates scat_sol */
623a5e57a09SHong Zhang       ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
624397b6df1SKris Buschelman     }
625a5e57a09SHong Zhang     if (!mumps->scat_sol) { /* create scatter scat_sol */
626a5e57a09SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,mumps->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */
627a5e57a09SHong Zhang       for (i=0; i<mumps->id.lsol_loc; i++) {
628a5e57a09SHong Zhang         mumps->id.isol_loc[i] -= 1; /* change Fortran style to C style */
629a5e57a09SHong Zhang       }
630a5e57a09SHong Zhang       ierr = ISCreateGeneral(PETSC_COMM_SELF,mumps->id.lsol_loc,mumps->id.isol_loc,PETSC_COPY_VALUES,&is_petsc);CHKERRQ(ierr);  /* to */
631a5e57a09SHong Zhang       ierr = VecScatterCreate(mumps->x_seq,is_iden,x,is_petsc,&mumps->scat_sol);CHKERRQ(ierr);
6326bf464f9SBarry Smith       ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
6336bf464f9SBarry Smith       ierr = ISDestroy(&is_petsc);CHKERRQ(ierr);
6342205254eSKarl Rupp 
635a5e57a09SHong Zhang       mumps->ICNTL9_pre = mumps->id.ICNTL(9); /* save current value of id.ICNTL(9) */
636397b6df1SKris Buschelman     }
637a5e57a09SHong Zhang 
638a5e57a09SHong Zhang     ierr = VecScatterBegin(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
639a5e57a09SHong Zhang     ierr = VecScatterEnd(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
640329ec9b3SHong Zhang   }
641397b6df1SKris Buschelman   PetscFunctionReturn(0);
642397b6df1SKris Buschelman }
643397b6df1SKris Buschelman 
64451d5961aSHong Zhang #undef __FUNCT__
64551d5961aSHong Zhang #define __FUNCT__ "MatSolveTranspose_MUMPS"
64651d5961aSHong Zhang PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x)
64751d5961aSHong Zhang {
648a5e57a09SHong Zhang   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
64951d5961aSHong Zhang   PetscErrorCode ierr;
65051d5961aSHong Zhang 
65151d5961aSHong Zhang   PetscFunctionBegin;
652a5e57a09SHong Zhang   mumps->id.ICNTL(9) = 0;
6530ad0caddSJed Brown   ierr = MatSolve_MUMPS(A,b,x);CHKERRQ(ierr);
654a5e57a09SHong Zhang   mumps->id.ICNTL(9) = 1;
65551d5961aSHong Zhang   PetscFunctionReturn(0);
65651d5961aSHong Zhang }
65751d5961aSHong Zhang 
658e0b74bf9SHong Zhang #undef __FUNCT__
659e0b74bf9SHong Zhang #define __FUNCT__ "MatMatSolve_MUMPS"
660e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X)
661e0b74bf9SHong Zhang {
662bda8bf91SBarry Smith   PetscErrorCode ierr;
663bda8bf91SBarry Smith   PetscBool      flg;
6644e34a73bSHong Zhang   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
665*71aed81dSHong Zhang   PetscInt       i,nrhs,m,M;
6662cd7d884SHong Zhang   PetscScalar    *array,*bray;
667bda8bf91SBarry Smith 
668e0b74bf9SHong Zhang   PetscFunctionBegin;
6690298fd71SBarry Smith   ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
670801fbe65SHong Zhang   if (!flg) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
6710298fd71SBarry Smith   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
672801fbe65SHong Zhang   if (!flg) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
673801fbe65SHong Zhang   if (B->rmap->n != X->rmap->n) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B and X must have same row distribution");
6744e34a73bSHong Zhang 
6752cd7d884SHong Zhang   ierr = MatGetLocalSize(B,&m,NULL);CHKERRQ(ierr);
6762cd7d884SHong Zhang   ierr = MatGetSize(B,&M,&nrhs);CHKERRQ(ierr);
6774e34a73bSHong Zhang 
6782cd7d884SHong Zhang   if (mumps->size == 1) {
6792cd7d884SHong Zhang     /* copy B to X */
6802cd7d884SHong Zhang     ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr);
6812cd7d884SHong Zhang     ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr);
6822cd7d884SHong Zhang     for (i=0; i<M*nrhs; i++) array[i] = bray[i];
6832cd7d884SHong Zhang     ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr);
6842cd7d884SHong Zhang 
6852cd7d884SHong Zhang     mumps->id.nrhs = nrhs;
6862cd7d884SHong Zhang     mumps->id.lrhs = M;
6872cd7d884SHong Zhang #if defined(PETSC_USE_COMPLEX)
6882cd7d884SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
6892cd7d884SHong Zhang     mumps->id.rhs = (mumps_complex*)array;
6902cd7d884SHong Zhang #else
6912cd7d884SHong Zhang     mumps->id.rhs = (mumps_double_complex*)array;
6922cd7d884SHong Zhang #endif
6932cd7d884SHong Zhang #else
6942cd7d884SHong Zhang     mumps->id.rhs = array;
6952cd7d884SHong Zhang #endif
696801fbe65SHong Zhang 
6972cd7d884SHong Zhang     /* solve phase */
6982cd7d884SHong Zhang     /*-------------*/
6992cd7d884SHong Zhang     mumps->id.job = JOB_SOLVE;
7002cd7d884SHong Zhang     PetscMUMPS_c(&mumps->id);
7012cd7d884SHong 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));
7022cd7d884SHong Zhang 
7032cd7d884SHong Zhang     ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr);
704801fbe65SHong Zhang   } else {  /************** parallel case ***************/
705*71aed81dSHong Zhang     PetscInt       lsol_loc,nlsol_loc,*isol_loc,*idx,*iidx,*idxx,*isol_loc_save;
706*71aed81dSHong Zhang     PetscScalar    *sol_loc,*sol_loc_save;
707801fbe65SHong Zhang     IS             is_to,is_from;
708801fbe65SHong Zhang     PetscInt       k,proc,j;
709801fbe65SHong Zhang     const PetscInt *rstart;
710*71aed81dSHong Zhang     Vec            v_mpi,bb_seq,x_seq;
71174f0fcc7SHong Zhang     VecScatter     scat_rhss, scat_sols;
712801fbe65SHong Zhang 
713801fbe65SHong Zhang     /* create x_seq to hold local solution */
714*71aed81dSHong Zhang     isol_loc_save = mumps->id.isol_loc; /* save it for MatSovle() */
715*71aed81dSHong Zhang     sol_loc_save  = mumps->id.sol_loc;
716801fbe65SHong Zhang 
717*71aed81dSHong Zhang     lsol_loc  = mumps->id.INFO(23);
718*71aed81dSHong Zhang     nlsol_loc = nrhs*lsol_loc;     /* length of sol_loc */
719*71aed81dSHong Zhang     ierr = PetscMalloc2(nlsol_loc,&sol_loc,nlsol_loc,&isol_loc);CHKERRQ(ierr);
720801fbe65SHong Zhang #if defined(PETSC_USE_COMPLEX)
721801fbe65SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
722801fbe65SHong Zhang     mumps->id.sol_loc = (mumps_complex*)sol_loc;
723801fbe65SHong Zhang #else
724801fbe65SHong Zhang     mumps->id.sol_loc = (mumps_double_complex*)sol_loc;
725801fbe65SHong Zhang #endif
726801fbe65SHong Zhang #else
727801fbe65SHong Zhang     mumps->id.sol_loc = sol_loc;
728801fbe65SHong Zhang #endif
729801fbe65SHong Zhang     mumps->id.isol_loc = isol_loc;
730801fbe65SHong Zhang 
731*71aed81dSHong Zhang     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,nlsol_loc,sol_loc,&x_seq);CHKERRQ(ierr);
7322cd7d884SHong Zhang 
73374f0fcc7SHong Zhang     /* copy rhs matrix B into vector v_mpi */
734801fbe65SHong Zhang     ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr);
73574f0fcc7SHong Zhang     ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)B),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi);CHKERRQ(ierr);
736801fbe65SHong Zhang     ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr);
737801fbe65SHong Zhang 
73874f0fcc7SHong Zhang     /* scatter v_mpi to bb_seq because MUMPS only supports centralized rhs */
73974f0fcc7SHong Zhang     /* idx: maps from k-th index of v_mpi to (i,j)-th global entry of B;
740801fbe65SHong Zhang       iidx: inverse of idx, will be used by scattering xx_seq -> X       */
741801fbe65SHong Zhang     ierr = PetscMalloc2(nrhs*M,&idx,nrhs*M,&iidx);CHKERRQ(ierr);
742801fbe65SHong Zhang     ierr = MatGetOwnershipRanges(B,&rstart);CHKERRQ(ierr);
743801fbe65SHong Zhang     k = 0;
744801fbe65SHong Zhang     for (proc=0; proc<mumps->size; proc++){
745801fbe65SHong Zhang       for (j=0; j<nrhs; j++){
746801fbe65SHong Zhang         for (i=rstart[proc]; i<rstart[proc+1]; i++){
747801fbe65SHong Zhang           iidx[j*M + i] = k;
748801fbe65SHong Zhang           idx[k++]      = j*M + i;
749801fbe65SHong Zhang         }
750801fbe65SHong Zhang       }
7512cd7d884SHong Zhang     }
7522cd7d884SHong Zhang 
753801fbe65SHong Zhang     if (!mumps->myid) {
75474f0fcc7SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,nrhs*M,&bb_seq);CHKERRQ(ierr);
755801fbe65SHong Zhang       ierr = ISCreateGeneral(PETSC_COMM_SELF,nrhs*M,idx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr);
756801fbe65SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,nrhs*M,0,1,&is_from);CHKERRQ(ierr);
757801fbe65SHong Zhang     } else {
75874f0fcc7SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&bb_seq);CHKERRQ(ierr);
759801fbe65SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_to);CHKERRQ(ierr);
760801fbe65SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_from);CHKERRQ(ierr);
761801fbe65SHong Zhang     }
76274f0fcc7SHong Zhang     ierr = VecScatterCreate(v_mpi,is_from,bb_seq,is_to,&scat_rhss);CHKERRQ(ierr);
76374f0fcc7SHong Zhang     ierr = VecScatterBegin(scat_rhss,v_mpi,bb_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
764801fbe65SHong Zhang     ierr = ISDestroy(&is_to);CHKERRQ(ierr);
765801fbe65SHong Zhang     ierr = ISDestroy(&is_from);CHKERRQ(ierr);
76674f0fcc7SHong Zhang     ierr = VecScatterEnd(scat_rhss,v_mpi,bb_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
767801fbe65SHong Zhang 
768801fbe65SHong Zhang     if (!mumps->myid) { /* define rhs on the host */
76974f0fcc7SHong Zhang       ierr = VecGetArray(bb_seq,&bray);CHKERRQ(ierr);
770801fbe65SHong Zhang       mumps->id.nrhs = nrhs;
771801fbe65SHong Zhang       mumps->id.lrhs = M;
772801fbe65SHong Zhang #if defined(PETSC_USE_COMPLEX)
773801fbe65SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
774801fbe65SHong Zhang       mumps->id.rhs = (mumps_complex*)bray;
775801fbe65SHong Zhang #else
776801fbe65SHong Zhang       mumps->id.rhs = (mumps_double_complex*)bray;
777801fbe65SHong Zhang #endif
778801fbe65SHong Zhang #else
779801fbe65SHong Zhang     mumps->id.rhs = bray;
780801fbe65SHong Zhang #endif
78174f0fcc7SHong Zhang     ierr = VecRestoreArray(bb_seq,&bray);CHKERRQ(ierr);
782801fbe65SHong Zhang     }
783801fbe65SHong Zhang 
784801fbe65SHong Zhang     /* solve phase */
785801fbe65SHong Zhang     /*-------------*/
786801fbe65SHong Zhang     mumps->id.job = JOB_SOLVE;
787801fbe65SHong Zhang     PetscMUMPS_c(&mumps->id);
788801fbe65SHong 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));
789801fbe65SHong Zhang 
790*71aed81dSHong Zhang     /* scatter mumps distributed solution to petsc vector xx_mpi, which shares local arrays with solution matrix X */
79174f0fcc7SHong Zhang     ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr);
79274f0fcc7SHong Zhang     ierr = VecPlaceArray(v_mpi,array);CHKERRQ(ierr);
793801fbe65SHong Zhang 
794801fbe65SHong Zhang     /* create scatter scat_sols */
795*71aed81dSHong Zhang     ierr = PetscMalloc1(nlsol_loc,&idxx);CHKERRQ(ierr);
796*71aed81dSHong Zhang     ierr = ISCreateStride(PETSC_COMM_SELF,nlsol_loc,0,1,&is_from);CHKERRQ(ierr);
797*71aed81dSHong Zhang     for (i=0; i<lsol_loc; i++) {
798801fbe65SHong Zhang       mumps->id.isol_loc[i] -= 1; /* change Fortran style to C style */
799801fbe65SHong Zhang       idxx[i] = iidx[mumps->id.isol_loc[i]];
800801fbe65SHong Zhang       for (j=1; j<nrhs; j++){
801*71aed81dSHong Zhang         idxx[j*lsol_loc+i] = iidx[mumps->id.isol_loc[i]+j*M];
802801fbe65SHong Zhang       }
803801fbe65SHong Zhang     }
804*71aed81dSHong Zhang     ierr = ISCreateGeneral(PETSC_COMM_SELF,nlsol_loc,idxx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr);
805*71aed81dSHong Zhang     ierr = VecScatterCreate(x_seq,is_from,v_mpi,is_to,&scat_sols);CHKERRQ(ierr);
806*71aed81dSHong Zhang     ierr = VecScatterBegin(scat_sols,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
807801fbe65SHong Zhang     ierr = ISDestroy(&is_from);CHKERRQ(ierr);
808801fbe65SHong Zhang     ierr = ISDestroy(&is_to);CHKERRQ(ierr);
809*71aed81dSHong Zhang     ierr = VecScatterEnd(scat_sols,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
810801fbe65SHong Zhang     ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr);
811*71aed81dSHong Zhang 
812*71aed81dSHong Zhang     /* free spaces */
813*71aed81dSHong Zhang     mumps->id.sol_loc = sol_loc_save;
814*71aed81dSHong Zhang     mumps->id.isol_loc = isol_loc_save;
815*71aed81dSHong Zhang 
816*71aed81dSHong Zhang     ierr = PetscFree2(sol_loc,isol_loc);CHKERRQ(ierr);
817801fbe65SHong Zhang     ierr = PetscFree2(idx,iidx);CHKERRQ(ierr);
818801fbe65SHong Zhang     ierr = PetscFree(idxx);CHKERRQ(ierr);
819*71aed81dSHong Zhang     ierr = VecDestroy(&x_seq);CHKERRQ(ierr);
82074f0fcc7SHong Zhang     ierr = VecDestroy(&v_mpi);CHKERRQ(ierr);
82174f0fcc7SHong Zhang     ierr = VecDestroy(&bb_seq);CHKERRQ(ierr);
82274f0fcc7SHong Zhang     ierr = VecScatterDestroy(&scat_rhss);CHKERRQ(ierr);
82374f0fcc7SHong Zhang     ierr = VecScatterDestroy(&scat_sols);CHKERRQ(ierr);
824801fbe65SHong Zhang   }
825e0b74bf9SHong Zhang   PetscFunctionReturn(0);
826e0b74bf9SHong Zhang }
827e0b74bf9SHong Zhang 
828ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX)
829a58c3f20SHong Zhang /*
830a58c3f20SHong Zhang   input:
831a58c3f20SHong Zhang    F:        numeric factor
832a58c3f20SHong Zhang   output:
833a58c3f20SHong Zhang    nneg:     total number of negative pivots
834a58c3f20SHong Zhang    nzero:    0
835a58c3f20SHong Zhang    npos:     (global dimension of F) - nneg
836a58c3f20SHong Zhang */
837a58c3f20SHong Zhang 
838a58c3f20SHong Zhang #undef __FUNCT__
839a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
840dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
841a58c3f20SHong Zhang {
842a5e57a09SHong Zhang   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
843dfbe8321SBarry Smith   PetscErrorCode ierr;
844c1490034SHong Zhang   PetscMPIInt    size;
845a58c3f20SHong Zhang 
846a58c3f20SHong Zhang   PetscFunctionBegin;
847ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)F),&size);CHKERRQ(ierr);
848bcb30aebSHong 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 */
849a5e57a09SHong 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));
850ed85ac9fSHong Zhang 
851710ac8efSHong Zhang   if (nneg) *nneg = mumps->id.INFOG(12);
852ed85ac9fSHong Zhang   if (nzero || npos) {
853ed85ac9fSHong 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");
854710ac8efSHong Zhang     if (nzero) *nzero = mumps->id.INFOG(28);
855710ac8efSHong Zhang     if (npos) *npos   = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28));
856a58c3f20SHong Zhang   }
857a58c3f20SHong Zhang   PetscFunctionReturn(0);
858a58c3f20SHong Zhang }
859ace3df97SHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */
860a58c3f20SHong Zhang 
861397b6df1SKris Buschelman #undef __FUNCT__
862f6c57405SHong Zhang #define __FUNCT__ "MatFactorNumeric_MUMPS"
8630481f469SBarry Smith PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
864af281ebdSHong Zhang {
865a5e57a09SHong Zhang   Mat_MUMPS      *mumps =(Mat_MUMPS*)(F)->spptr;
8666849ba73SBarry Smith   PetscErrorCode ierr;
867e09efc27SHong Zhang   Mat            F_diag;
868ace3abfcSBarry Smith   PetscBool      isMPIAIJ;
869397b6df1SKris Buschelman 
870397b6df1SKris Buschelman   PetscFunctionBegin;
871a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
872397b6df1SKris Buschelman 
873397b6df1SKris Buschelman   /* numerical factorization phase */
874329ec9b3SHong Zhang   /*-------------------------------*/
875a5e57a09SHong Zhang   mumps->id.job = JOB_FACTNUMERIC;
8764e34a73bSHong Zhang   if (!mumps->id.ICNTL(18)) { /* A is centralized */
877a5e57a09SHong Zhang     if (!mumps->myid) {
878397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
8792907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
880a5e57a09SHong Zhang       mumps->id.a = (mumps_complex*)mumps->val;
8812907cef9SHong Zhang #else
882a5e57a09SHong Zhang       mumps->id.a = (mumps_double_complex*)mumps->val;
8832907cef9SHong Zhang #endif
884397b6df1SKris Buschelman #else
885a5e57a09SHong Zhang       mumps->id.a = mumps->val;
886397b6df1SKris Buschelman #endif
887397b6df1SKris Buschelman     }
888397b6df1SKris Buschelman   } else {
889397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
8902907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
891a5e57a09SHong Zhang     mumps->id.a_loc = (mumps_complex*)mumps->val;
8922907cef9SHong Zhang #else
893a5e57a09SHong Zhang     mumps->id.a_loc = (mumps_double_complex*)mumps->val;
8942907cef9SHong Zhang #endif
895397b6df1SKris Buschelman #else
896a5e57a09SHong Zhang     mumps->id.a_loc = mumps->val;
897397b6df1SKris Buschelman #endif
898397b6df1SKris Buschelman   }
899a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
900a5e57a09SHong Zhang   if (mumps->id.INFOG(1) < 0) {
901151787a6SHong Zhang     if (mumps->id.INFO(1) == -13) {
902151787a6SHong Zhang       if (mumps->id.INFO(2) < 0) {
903151787a6SHong 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));
904151787a6SHong Zhang       } else {
905151787a6SHong 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));
906151787a6SHong Zhang       }
907151787a6SHong 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));
908397b6df1SKris Buschelman   }
909a5e57a09SHong 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));
910397b6df1SKris Buschelman 
911dcd589f8SShri Abhyankar   (F)->assembled      = PETSC_TRUE;
912a5e57a09SHong Zhang   mumps->matstruc     = SAME_NONZERO_PATTERN;
913a5e57a09SHong Zhang   mumps->CleanUpMUMPS = PETSC_TRUE;
91467877ebaSShri Abhyankar 
915a5e57a09SHong Zhang   if (mumps->size > 1) {
91667877ebaSShri Abhyankar     PetscInt    lsol_loc;
91767877ebaSShri Abhyankar     PetscScalar *sol_loc;
9182205254eSKarl Rupp 
919c2093ab7SHong Zhang     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr);
920c2093ab7SHong Zhang     if (isMPIAIJ) F_diag = ((Mat_MPIAIJ*)(F)->data)->A;
921c2093ab7SHong Zhang     else F_diag = ((Mat_MPISBAIJ*)(F)->data)->A;
922c2093ab7SHong Zhang     F_diag->assembled = PETSC_TRUE;
923c2093ab7SHong Zhang 
924c2093ab7SHong Zhang     /* distributed solution; Create x_seq=sol_loc for repeated use */
925c2093ab7SHong Zhang     if (mumps->x_seq) {
926c2093ab7SHong Zhang       ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
927c2093ab7SHong Zhang       ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr);
928c2093ab7SHong Zhang       ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr);
929c2093ab7SHong Zhang     }
930a5e57a09SHong Zhang     lsol_loc = mumps->id.INFO(23); /* length of sol_loc */
931dcca6d9dSJed Brown     ierr = PetscMalloc2(lsol_loc,&sol_loc,lsol_loc,&mumps->id.isol_loc);CHKERRQ(ierr);
932a5e57a09SHong Zhang     mumps->id.lsol_loc = lsol_loc;
93367877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
9342907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
935a5e57a09SHong Zhang     mumps->id.sol_loc = (mumps_complex*)sol_loc;
9362907cef9SHong Zhang #else
937a5e57a09SHong Zhang     mumps->id.sol_loc = (mumps_double_complex*)sol_loc;
9382907cef9SHong Zhang #endif
93967877ebaSShri Abhyankar #else
940a5e57a09SHong Zhang     mumps->id.sol_loc = sol_loc;
94167877ebaSShri Abhyankar #endif
942a5e57a09SHong Zhang     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq);CHKERRQ(ierr);
94367877ebaSShri Abhyankar   }
944397b6df1SKris Buschelman   PetscFunctionReturn(0);
945397b6df1SKris Buschelman }
946397b6df1SKris Buschelman 
9479a2535b5SHong Zhang /* Sets MUMPS options from the options database */
948dcd589f8SShri Abhyankar #undef __FUNCT__
9499a2535b5SHong Zhang #define __FUNCT__ "PetscSetMUMPSFromOptions"
9509a2535b5SHong Zhang PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A)
951dcd589f8SShri Abhyankar {
9529a2535b5SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
953dcd589f8SShri Abhyankar   PetscErrorCode ierr;
954dcd589f8SShri Abhyankar   PetscInt       icntl;
955ace3abfcSBarry Smith   PetscBool      flg;
956dcd589f8SShri Abhyankar 
957dcd589f8SShri Abhyankar   PetscFunctionBegin;
958ce94432eSBarry Smith   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
9599a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);CHKERRQ(ierr);
9609a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(1) = icntl;
9619a2535b5SHong 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);
9629a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(2) = icntl;
9639a2535b5SHong 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);
9649a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(3) = icntl;
965dcd589f8SShri Abhyankar 
9669a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",mumps->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
9679a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(4) = icntl;
9689a2535b5SHong Zhang   if (mumps->id.ICNTL(4) || PetscLogPrintInfo) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */
9699a2535b5SHong Zhang 
970d341cd04SHong 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);
9719a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(6) = icntl;
9729a2535b5SHong Zhang 
973d341cd04SHong 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);
974dcd589f8SShri Abhyankar   if (flg) {
9752205254eSKarl 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");
9762205254eSKarl Rupp     else mumps->id.ICNTL(7) = icntl;
977dcd589f8SShri Abhyankar   }
978e0b74bf9SHong Zhang 
9790298fd71SBarry 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);
980d341cd04SHong 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() */
9810298fd71SBarry 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);
982d341cd04SHong 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);
983d341cd04SHong 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);
984d341cd04SHong 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);
985d341cd04SHong 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);
986d341cd04SHong 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);
9874e34a73bSHong 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 */
988d341cd04SHong 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 */
9899a2535b5SHong Zhang 
990d341cd04SHong 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);
9910298fd71SBarry 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);
9920298fd71SBarry 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);
9939a2535b5SHong Zhang   if (mumps->id.ICNTL(24)) {
9949a2535b5SHong Zhang     mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */
995d7ebd59bSHong Zhang   }
996d7ebd59bSHong Zhang 
997d341cd04SHong 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);
998d341cd04SHong 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);
9992cd7d884SHong 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);
10000298fd71SBarry 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);
1001d341cd04SHong 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);
10020298fd71SBarry 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);
1003d341cd04SHong 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);
10044e34a73bSHong 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 */
10050298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),NULL);CHKERRQ(ierr);
1006dcd589f8SShri Abhyankar 
10070298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),NULL);CHKERRQ(ierr);
10080298fd71SBarry 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);
10090298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),NULL);CHKERRQ(ierr);
10100298fd71SBarry 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);
10110298fd71SBarry 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);
1012e5bb22a1SHong Zhang 
10130298fd71SBarry Smith   ierr = PetscOptionsString("-mat_mumps_ooc_tmpdir", "out of core directory", "None", mumps->id.ooc_tmpdir, mumps->id.ooc_tmpdir, 256, NULL);
1014dcd589f8SShri Abhyankar   PetscOptionsEnd();
1015dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
1016dcd589f8SShri Abhyankar }
1017dcd589f8SShri Abhyankar 
1018dcd589f8SShri Abhyankar #undef __FUNCT__
1019dcd589f8SShri Abhyankar #define __FUNCT__ "PetscInitializeMUMPS"
1020f697e70eSHong Zhang PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps)
1021dcd589f8SShri Abhyankar {
1022dcd589f8SShri Abhyankar   PetscErrorCode ierr;
1023dcd589f8SShri Abhyankar 
1024dcd589f8SShri Abhyankar   PetscFunctionBegin;
1025ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)A), &mumps->myid);
1026ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&mumps->size);CHKERRQ(ierr);
1027ce94432eSBarry Smith   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(mumps->comm_mumps));CHKERRQ(ierr);
10282205254eSKarl Rupp 
1029f697e70eSHong Zhang   mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps);
1030f697e70eSHong Zhang 
1031f697e70eSHong Zhang   mumps->id.job = JOB_INIT;
1032f697e70eSHong Zhang   mumps->id.par = 1;  /* host participates factorizaton and solve */
1033f697e70eSHong Zhang   mumps->id.sym = mumps->sym;
10342907cef9SHong Zhang   PetscMUMPS_c(&mumps->id);
1035f697e70eSHong Zhang 
1036f697e70eSHong Zhang   mumps->CleanUpMUMPS = PETSC_FALSE;
10370298fd71SBarry Smith   mumps->scat_rhs     = NULL;
10380298fd71SBarry Smith   mumps->scat_sol     = NULL;
10399a2535b5SHong Zhang 
104070544d5fSHong Zhang   /* set PETSc-MUMPS default options - override MUMPS default */
10419a2535b5SHong Zhang   mumps->id.ICNTL(3) = 0;
10429a2535b5SHong Zhang   mumps->id.ICNTL(4) = 0;
10439a2535b5SHong Zhang   if (mumps->size == 1) {
10449a2535b5SHong Zhang     mumps->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
10459a2535b5SHong Zhang   } else {
10469a2535b5SHong Zhang     mumps->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
10474e34a73bSHong Zhang     mumps->id.ICNTL(20) = 0;   /* rhs is in dense format */
104870544d5fSHong Zhang     mumps->id.ICNTL(21) = 1;   /* distributed solution */
10499a2535b5SHong Zhang   }
1050dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
1051dcd589f8SShri Abhyankar }
1052dcd589f8SShri Abhyankar 
1053a5e57a09SHong 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 */
1054397b6df1SKris Buschelman #undef __FUNCT__
1055f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
10560481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1057b24902e0SBarry Smith {
1058a5e57a09SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
1059dcd589f8SShri Abhyankar   PetscErrorCode ierr;
106067877ebaSShri Abhyankar   Vec            b;
106167877ebaSShri Abhyankar   IS             is_iden;
106267877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1063397b6df1SKris Buschelman 
1064397b6df1SKris Buschelman   PetscFunctionBegin;
1065a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1066dcd589f8SShri Abhyankar 
10679a2535b5SHong Zhang   /* Set MUMPS options from the options database */
10689a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1069dcd589f8SShri Abhyankar 
1070a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1071dcd589f8SShri Abhyankar 
107267877ebaSShri Abhyankar   /* analysis phase */
107367877ebaSShri Abhyankar   /*----------------*/
1074a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1075a5e57a09SHong Zhang   mumps->id.n   = M;
1076a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
107767877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1078a5e57a09SHong Zhang     if (!mumps->myid) {
1079a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1080a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
108167877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
10822907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1083a5e57a09SHong Zhang         mumps->id.a = (mumps_complex*)mumps->val;
10842907cef9SHong Zhang #else
1085a5e57a09SHong Zhang         mumps->id.a = (mumps_double_complex*)mumps->val;
10862907cef9SHong Zhang #endif
108767877ebaSShri Abhyankar #else
1088a5e57a09SHong Zhang         mumps->id.a = mumps->val;
108967877ebaSShri Abhyankar #endif
109067877ebaSShri Abhyankar       }
1091a5e57a09SHong Zhang       if (mumps->id.ICNTL(7) == 1) { /* use user-provide matrix ordering - assuming r = c ordering */
10925248a706SHong Zhang         /*
10935248a706SHong Zhang         PetscBool      flag;
10945248a706SHong Zhang         ierr = ISEqual(r,c,&flag);CHKERRQ(ierr);
10955248a706SHong Zhang         if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"row_perm != col_perm");
10965248a706SHong Zhang         ierr = ISView(r,PETSC_VIEWER_STDOUT_SELF);
10975248a706SHong Zhang          */
1098a5e57a09SHong Zhang         if (!mumps->myid) {
1099e0b74bf9SHong Zhang           const PetscInt *idx;
1100e0b74bf9SHong Zhang           PetscInt       i,*perm_in;
11012205254eSKarl Rupp 
1102785e854fSJed Brown           ierr = PetscMalloc1(M,&perm_in);CHKERRQ(ierr);
1103e0b74bf9SHong Zhang           ierr = ISGetIndices(r,&idx);CHKERRQ(ierr);
11042205254eSKarl Rupp 
1105a5e57a09SHong Zhang           mumps->id.perm_in = perm_in;
1106e0b74bf9SHong Zhang           for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */
1107e0b74bf9SHong Zhang           ierr = ISRestoreIndices(r,&idx);CHKERRQ(ierr);
1108e0b74bf9SHong Zhang         }
1109e0b74bf9SHong Zhang       }
111067877ebaSShri Abhyankar     }
111167877ebaSShri Abhyankar     break;
111267877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1113a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1114a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1115a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
111667877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
11172907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1118a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_complex*)mumps->val;
11192907cef9SHong Zhang #else
1120a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_double_complex*)mumps->val;
11212907cef9SHong Zhang #endif
112267877ebaSShri Abhyankar #else
1123a5e57a09SHong Zhang       mumps->id.a_loc = mumps->val;
112467877ebaSShri Abhyankar #endif
112567877ebaSShri Abhyankar     }
112667877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1127a5e57a09SHong Zhang     if (!mumps->myid) {
11282cd7d884SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->rmap->N,&mumps->b_seq);CHKERRQ(ierr);
11292cd7d884SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->N,0,1,&is_iden);CHKERRQ(ierr);
113067877ebaSShri Abhyankar     } else {
1131a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
113267877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
113367877ebaSShri Abhyankar     }
11342a7a6963SBarry Smith     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1135a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
11366bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
11376bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
113867877ebaSShri Abhyankar     break;
113967877ebaSShri Abhyankar   }
1140a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
1141a5e57a09SHong 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));
114267877ebaSShri Abhyankar 
1143719d5645SBarry Smith   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1144dcd589f8SShri Abhyankar   F->ops->solve           = MatSolve_MUMPS;
114551d5961aSHong Zhang   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
11464e34a73bSHong Zhang   F->ops->matsolve        = MatMatSolve_MUMPS;
1147b24902e0SBarry Smith   PetscFunctionReturn(0);
1148b24902e0SBarry Smith }
1149b24902e0SBarry Smith 
1150450b117fSShri Abhyankar /* Note the Petsc r and c permutations are ignored */
1151450b117fSShri Abhyankar #undef __FUNCT__
1152450b117fSShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_BAIJMUMPS"
1153450b117fSShri Abhyankar PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1154450b117fSShri Abhyankar {
1155a5e57a09SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
1156dcd589f8SShri Abhyankar   PetscErrorCode ierr;
115767877ebaSShri Abhyankar   Vec            b;
115867877ebaSShri Abhyankar   IS             is_iden;
115967877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1160450b117fSShri Abhyankar 
1161450b117fSShri Abhyankar   PetscFunctionBegin;
1162a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1163dcd589f8SShri Abhyankar 
11649a2535b5SHong Zhang   /* Set MUMPS options from the options database */
11659a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1166dcd589f8SShri Abhyankar 
1167a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
116867877ebaSShri Abhyankar 
116967877ebaSShri Abhyankar   /* analysis phase */
117067877ebaSShri Abhyankar   /*----------------*/
1171a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1172a5e57a09SHong Zhang   mumps->id.n   = M;
1173a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
117467877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1175a5e57a09SHong Zhang     if (!mumps->myid) {
1176a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=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 = (mumps_complex*)mumps->val;
11812907cef9SHong Zhang #else
1182a5e57a09SHong Zhang         mumps->id.a = (mumps_double_complex*)mumps->val;
11832907cef9SHong Zhang #endif
118467877ebaSShri Abhyankar #else
1185a5e57a09SHong Zhang         mumps->id.a = mumps->val;
118667877ebaSShri Abhyankar #endif
118767877ebaSShri Abhyankar       }
118867877ebaSShri Abhyankar     }
118967877ebaSShri Abhyankar     break;
119067877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1191a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1192a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1193a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
119467877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
11952907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1196a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_complex*)mumps->val;
11972907cef9SHong Zhang #else
1198a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_double_complex*)mumps->val;
11992907cef9SHong Zhang #endif
120067877ebaSShri Abhyankar #else
1201a5e57a09SHong Zhang       mumps->id.a_loc = mumps->val;
120267877ebaSShri Abhyankar #endif
120367877ebaSShri Abhyankar     }
120467877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1205a5e57a09SHong Zhang     if (!mumps->myid) {
1206a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
120767877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
120867877ebaSShri Abhyankar     } else {
1209a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
121067877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
121167877ebaSShri Abhyankar     }
12122a7a6963SBarry Smith     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1213a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
12146bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
12156bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
121667877ebaSShri Abhyankar     break;
121767877ebaSShri Abhyankar   }
1218a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
1219a5e57a09SHong 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));
122067877ebaSShri Abhyankar 
1221450b117fSShri Abhyankar   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1222dcd589f8SShri Abhyankar   F->ops->solve           = MatSolve_MUMPS;
122351d5961aSHong Zhang   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1224450b117fSShri Abhyankar   PetscFunctionReturn(0);
1225450b117fSShri Abhyankar }
1226b24902e0SBarry Smith 
1227141f4205SHong Zhang /* Note the Petsc r permutation and factor info are ignored */
1228397b6df1SKris Buschelman #undef __FUNCT__
122967877ebaSShri Abhyankar #define __FUNCT__ "MatCholeskyFactorSymbolic_MUMPS"
123067877ebaSShri Abhyankar PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
1231b24902e0SBarry Smith {
1232a5e57a09SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
1233dcd589f8SShri Abhyankar   PetscErrorCode ierr;
123467877ebaSShri Abhyankar   Vec            b;
123567877ebaSShri Abhyankar   IS             is_iden;
123667877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1237397b6df1SKris Buschelman 
1238397b6df1SKris Buschelman   PetscFunctionBegin;
1239a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1240dcd589f8SShri Abhyankar 
12419a2535b5SHong Zhang   /* Set MUMPS options from the options database */
12429a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1243dcd589f8SShri Abhyankar 
1244a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1245dcd589f8SShri Abhyankar 
124667877ebaSShri Abhyankar   /* analysis phase */
124767877ebaSShri Abhyankar   /*----------------*/
1248a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1249a5e57a09SHong Zhang   mumps->id.n   = M;
1250a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
125167877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1252a5e57a09SHong Zhang     if (!mumps->myid) {
1253a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1254a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
125567877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
12562907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1257a5e57a09SHong Zhang         mumps->id.a = (mumps_complex*)mumps->val;
12582907cef9SHong Zhang #else
1259a5e57a09SHong Zhang         mumps->id.a = (mumps_double_complex*)mumps->val;
12602907cef9SHong Zhang #endif
126167877ebaSShri Abhyankar #else
1262a5e57a09SHong Zhang         mumps->id.a = mumps->val;
126367877ebaSShri Abhyankar #endif
126467877ebaSShri Abhyankar       }
126567877ebaSShri Abhyankar     }
126667877ebaSShri Abhyankar     break;
126767877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1268a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1269a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1270a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
127167877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
12722907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1273a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_complex*)mumps->val;
12742907cef9SHong Zhang #else
1275a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_double_complex*)mumps->val;
12762907cef9SHong Zhang #endif
127767877ebaSShri Abhyankar #else
1278a5e57a09SHong Zhang       mumps->id.a_loc = mumps->val;
127967877ebaSShri Abhyankar #endif
128067877ebaSShri Abhyankar     }
128167877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1282a5e57a09SHong Zhang     if (!mumps->myid) {
1283a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
128467877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
128567877ebaSShri Abhyankar     } else {
1286a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
128767877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
128867877ebaSShri Abhyankar     }
12892a7a6963SBarry Smith     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1290a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
12916bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
12926bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
129367877ebaSShri Abhyankar     break;
129467877ebaSShri Abhyankar   }
1295a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
1296a5e57a09SHong 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));
129767877ebaSShri Abhyankar 
12982792810eSHong Zhang   F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
1299dcd589f8SShri Abhyankar   F->ops->solve                 = MatSolve_MUMPS;
130051d5961aSHong Zhang   F->ops->solvetranspose        = MatSolve_MUMPS;
13014e34a73bSHong Zhang   F->ops->matsolve              = MatMatSolve_MUMPS;
13024e34a73bSHong Zhang #if defined(PETSC_USE_COMPLEX)
13030298fd71SBarry Smith   F->ops->getinertia = NULL;
13044e34a73bSHong Zhang #else
13054e34a73bSHong Zhang   F->ops->getinertia = MatGetInertia_SBAIJMUMPS;
1306db4efbfdSBarry Smith #endif
1307b24902e0SBarry Smith   PetscFunctionReturn(0);
1308b24902e0SBarry Smith }
1309b24902e0SBarry Smith 
13104e34a73bSHong Zhang //update!!!
1311397b6df1SKris Buschelman #undef __FUNCT__
131264e6c443SBarry Smith #define __FUNCT__ "MatView_MUMPS"
131364e6c443SBarry Smith PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
131474ed9c26SBarry Smith {
1315f6c57405SHong Zhang   PetscErrorCode    ierr;
131664e6c443SBarry Smith   PetscBool         iascii;
131764e6c443SBarry Smith   PetscViewerFormat format;
1318a5e57a09SHong Zhang   Mat_MUMPS         *mumps=(Mat_MUMPS*)A->spptr;
1319f6c57405SHong Zhang 
1320f6c57405SHong Zhang   PetscFunctionBegin;
132164e6c443SBarry Smith   /* check if matrix is mumps type */
132264e6c443SBarry Smith   if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0);
132364e6c443SBarry Smith 
1324251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
132564e6c443SBarry Smith   if (iascii) {
132664e6c443SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
132764e6c443SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO) {
132864e6c443SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
1329a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                   %d \n",mumps->id.sym);CHKERRQ(ierr);
1330a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):            %d \n",mumps->id.par);CHKERRQ(ierr);
1331a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):         %d \n",mumps->id.ICNTL(1));CHKERRQ(ierr);
1332a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg): %d \n",mumps->id.ICNTL(2));CHKERRQ(ierr);
1333a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):   %d \n",mumps->id.ICNTL(3));CHKERRQ(ierr);
1334a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):        %d \n",mumps->id.ICNTL(4));CHKERRQ(ierr);
1335a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):         %d \n",mumps->id.ICNTL(5));CHKERRQ(ierr);
1336a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):        %d \n",mumps->id.ICNTL(6));CHKERRQ(ierr);
1337a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (sequentia matrix ordering):%d \n",mumps->id.ICNTL(7));CHKERRQ(ierr);
1338a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scalling strategy):        %d \n",mumps->id.ICNTL(8));CHKERRQ(ierr);
1339a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements):  %d \n",mumps->id.ICNTL(10));CHKERRQ(ierr);
1340a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):          %d \n",mumps->id.ICNTL(11));CHKERRQ(ierr);
1341a5e57a09SHong Zhang       if (mumps->id.ICNTL(11)>0) {
1342a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(4) (inf norm of input mat):        %g\n",mumps->id.RINFOG(4));CHKERRQ(ierr);
1343a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(5) (inf norm of solution):         %g\n",mumps->id.RINFOG(5));CHKERRQ(ierr);
1344a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(6) (inf norm of residual):         %g\n",mumps->id.RINFOG(6));CHKERRQ(ierr);
1345a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8));CHKERRQ(ierr);
1346a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(9) (error estimate):               %g \n",mumps->id.RINFOG(9));CHKERRQ(ierr);
1347a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11));CHKERRQ(ierr);
1348f6c57405SHong Zhang       }
1349a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",mumps->id.ICNTL(12));CHKERRQ(ierr);
1350a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",mumps->id.ICNTL(13));CHKERRQ(ierr);
1351a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",mumps->id.ICNTL(14));CHKERRQ(ierr);
1352f6c57405SHong Zhang       /* ICNTL(15-17) not used */
1353a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",mumps->id.ICNTL(18));CHKERRQ(ierr);
1354a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Shur complement info):                       %d \n",mumps->id.ICNTL(19));CHKERRQ(ierr);
1355a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",mumps->id.ICNTL(20));CHKERRQ(ierr);
1356a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (somumpstion struct):                            %d \n",mumps->id.ICNTL(21));CHKERRQ(ierr);
1357a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",mumps->id.ICNTL(22));CHKERRQ(ierr);
1358a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",mumps->id.ICNTL(23));CHKERRQ(ierr);
1359c0165424SHong Zhang 
1360a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",mumps->id.ICNTL(24));CHKERRQ(ierr);
1361a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",mumps->id.ICNTL(25));CHKERRQ(ierr);
1362a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",mumps->id.ICNTL(26));CHKERRQ(ierr);
1363a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",mumps->id.ICNTL(27));CHKERRQ(ierr);
1364a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(28) (use parallel or sequential ordering):        %d \n",mumps->id.ICNTL(28));CHKERRQ(ierr);
1365a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(29) (parallel ordering):                          %d \n",mumps->id.ICNTL(29));CHKERRQ(ierr);
136642179a6aSHong Zhang 
1367a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(30) (user-specified set of entries in inv(A)):    %d \n",mumps->id.ICNTL(30));CHKERRQ(ierr);
1368a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(31) (factors is discarded in the solve phase):    %d \n",mumps->id.ICNTL(31));CHKERRQ(ierr);
1369a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(33) (compute determinant):                        %d \n",mumps->id.ICNTL(33));CHKERRQ(ierr);
1370f6c57405SHong Zhang 
1371a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",mumps->id.CNTL(1));CHKERRQ(ierr);
1372a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2));CHKERRQ(ierr);
1373a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absomumpste pivoting threshold):      %g \n",mumps->id.CNTL(3));CHKERRQ(ierr);
1374a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(4) (vamumpse of static pivoting):         %g \n",mumps->id.CNTL(4));CHKERRQ(ierr);
1375a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",mumps->id.CNTL(5));CHKERRQ(ierr);
1376f6c57405SHong Zhang 
1377f6c57405SHong Zhang       /* infomation local to each processor */
137834ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);
13797b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
1380a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %g \n",mumps->myid,mumps->id.RINFO(1));CHKERRQ(ierr);
138134ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
138234ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);
1383a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(2));CHKERRQ(ierr);
138434ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
138534ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);
1386a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(3));CHKERRQ(ierr);
138734ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1388f6c57405SHong Zhang 
138934ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr);
1390a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] %d \n",mumps->myid,mumps->id.INFO(15));CHKERRQ(ierr);
139134ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1392f6c57405SHong Zhang 
139334ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr);
1394a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(16));CHKERRQ(ierr);
139534ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1396f6c57405SHong Zhang 
139734ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);
1398a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(23));CHKERRQ(ierr);
139934ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
14007b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
1401f6c57405SHong Zhang 
1402a5e57a09SHong Zhang       if (!mumps->myid) { /* information from the host */
1403a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",mumps->id.RINFOG(1));CHKERRQ(ierr);
1404a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",mumps->id.RINFOG(2));CHKERRQ(ierr);
1405a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",mumps->id.RINFOG(3));CHKERRQ(ierr);
1406a5e57a09SHong 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);
1407f6c57405SHong Zhang 
1408a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(3));CHKERRQ(ierr);
1409a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(4));CHKERRQ(ierr);
1410a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",mumps->id.INFOG(5));CHKERRQ(ierr);
1411a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",mumps->id.INFOG(6));CHKERRQ(ierr);
1412a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively use after analysis): %d \n",mumps->id.INFOG(7));CHKERRQ(ierr);
1413a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",mumps->id.INFOG(8));CHKERRQ(ierr);
1414a5e57a09SHong 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);
1415a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",mumps->id.INFOG(10));CHKERRQ(ierr);
1416a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",mumps->id.INFOG(11));CHKERRQ(ierr);
1417a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",mumps->id.INFOG(12));CHKERRQ(ierr);
1418a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",mumps->id.INFOG(13));CHKERRQ(ierr);
1419a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",mumps->id.INFOG(14));CHKERRQ(ierr);
1420a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",mumps->id.INFOG(15));CHKERRQ(ierr);
1421a5e57a09SHong 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);
1422a5e57a09SHong 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);
1423a5e57a09SHong 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);
1424a5e57a09SHong 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);
1425a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",mumps->id.INFOG(20));CHKERRQ(ierr);
1426a5e57a09SHong 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);
1427a5e57a09SHong 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);
1428a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",mumps->id.INFOG(23));CHKERRQ(ierr);
1429a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",mumps->id.INFOG(24));CHKERRQ(ierr);
1430a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",mumps->id.INFOG(25));CHKERRQ(ierr);
143140d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(28) (after factorization: number of null pivots encountered): %d\n",mumps->id.INFOG(28));CHKERRQ(ierr);
143240d435e3SHong 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);
143340d435e3SHong 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);
143440d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(32) (after analysis: type of analysis done): %d\n",mumps->id.INFOG(32));CHKERRQ(ierr);
143540d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(33) (value used for ICNTL(8)): %d\n",mumps->id.INFOG(33));CHKERRQ(ierr);
143640d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(34) (exponent of the determinant if determinant is requested): %d\n",mumps->id.INFOG(34));CHKERRQ(ierr);
1437f6c57405SHong Zhang       }
1438f6c57405SHong Zhang     }
1439cb828f0fSHong Zhang   }
1440f6c57405SHong Zhang   PetscFunctionReturn(0);
1441f6c57405SHong Zhang }
1442f6c57405SHong Zhang 
144335bd34faSBarry Smith #undef __FUNCT__
144435bd34faSBarry Smith #define __FUNCT__ "MatGetInfo_MUMPS"
144535bd34faSBarry Smith PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
144635bd34faSBarry Smith {
1447cb828f0fSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)A->spptr;
144835bd34faSBarry Smith 
144935bd34faSBarry Smith   PetscFunctionBegin;
145035bd34faSBarry Smith   info->block_size        = 1.0;
1451cb828f0fSHong Zhang   info->nz_allocated      = mumps->id.INFOG(20);
1452cb828f0fSHong Zhang   info->nz_used           = mumps->id.INFOG(20);
145335bd34faSBarry Smith   info->nz_unneeded       = 0.0;
145435bd34faSBarry Smith   info->assemblies        = 0.0;
145535bd34faSBarry Smith   info->mallocs           = 0.0;
145635bd34faSBarry Smith   info->memory            = 0.0;
145735bd34faSBarry Smith   info->fill_ratio_given  = 0;
145835bd34faSBarry Smith   info->fill_ratio_needed = 0;
145935bd34faSBarry Smith   info->factor_mallocs    = 0;
146035bd34faSBarry Smith   PetscFunctionReturn(0);
146135bd34faSBarry Smith }
146235bd34faSBarry Smith 
14635ccb76cbSHong Zhang /* -------------------------------------------------------------------------------------------*/
14645ccb76cbSHong Zhang #undef __FUNCT__
14655ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl_MUMPS"
14665ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival)
14675ccb76cbSHong Zhang {
1468a5e57a09SHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
14695ccb76cbSHong Zhang 
14705ccb76cbSHong Zhang   PetscFunctionBegin;
1471a5e57a09SHong Zhang   mumps->id.ICNTL(icntl) = ival;
14725ccb76cbSHong Zhang   PetscFunctionReturn(0);
14735ccb76cbSHong Zhang }
14745ccb76cbSHong Zhang 
14755ccb76cbSHong Zhang #undef __FUNCT__
1476bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetIcntl_MUMPS"
1477bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt *ival)
1478bc6112feSHong Zhang {
1479bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1480bc6112feSHong Zhang 
1481bc6112feSHong Zhang   PetscFunctionBegin;
1482bc6112feSHong Zhang   *ival = mumps->id.ICNTL(icntl);
1483bc6112feSHong Zhang   PetscFunctionReturn(0);
1484bc6112feSHong Zhang }
1485bc6112feSHong Zhang 
1486bc6112feSHong Zhang #undef __FUNCT__
14875ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl"
14885ccb76cbSHong Zhang /*@
14895ccb76cbSHong Zhang   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()
14905ccb76cbSHong Zhang 
14915ccb76cbSHong Zhang    Logically Collective on Mat
14925ccb76cbSHong Zhang 
14935ccb76cbSHong Zhang    Input Parameters:
14945ccb76cbSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
14955ccb76cbSHong Zhang .  icntl - index of MUMPS parameter array ICNTL()
14965ccb76cbSHong Zhang -  ival - value of MUMPS ICNTL(icntl)
14975ccb76cbSHong Zhang 
14985ccb76cbSHong Zhang   Options Database:
14995ccb76cbSHong Zhang .   -mat_mumps_icntl_<icntl> <ival>
15005ccb76cbSHong Zhang 
15015ccb76cbSHong Zhang    Level: beginner
15025ccb76cbSHong Zhang 
15035ccb76cbSHong Zhang    References: MUMPS Users' Guide
15045ccb76cbSHong Zhang 
15055ccb76cbSHong Zhang .seealso: MatGetFactor()
15065ccb76cbSHong Zhang @*/
15075ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival)
15085ccb76cbSHong Zhang {
15095ccb76cbSHong Zhang   PetscErrorCode ierr;
15105ccb76cbSHong Zhang 
15115ccb76cbSHong Zhang   PetscFunctionBegin;
15125ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
15135ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,ival,3);
15145ccb76cbSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr);
15155ccb76cbSHong Zhang   PetscFunctionReturn(0);
15165ccb76cbSHong Zhang }
15175ccb76cbSHong Zhang 
1518bc6112feSHong Zhang #undef __FUNCT__
1519bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetIcntl"
1520a21f80fcSHong Zhang /*@
1521a21f80fcSHong Zhang   MatMumpsGetIcntl - Get MUMPS parameter ICNTL()
1522a21f80fcSHong Zhang 
1523a21f80fcSHong Zhang    Logically Collective on Mat
1524a21f80fcSHong Zhang 
1525a21f80fcSHong Zhang    Input Parameters:
1526a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1527a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array ICNTL()
1528a21f80fcSHong Zhang 
1529a21f80fcSHong Zhang   Output Parameter:
1530a21f80fcSHong Zhang .  ival - value of MUMPS ICNTL(icntl)
1531a21f80fcSHong Zhang 
1532a21f80fcSHong Zhang    Level: beginner
1533a21f80fcSHong Zhang 
1534a21f80fcSHong Zhang    References: MUMPS Users' Guide
1535a21f80fcSHong Zhang 
1536a21f80fcSHong Zhang .seealso: MatGetFactor()
1537a21f80fcSHong Zhang @*/
1538bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl(Mat F,PetscInt icntl,PetscInt *ival)
1539bc6112feSHong Zhang {
1540bc6112feSHong Zhang   PetscErrorCode ierr;
1541bc6112feSHong Zhang 
1542bc6112feSHong Zhang   PetscFunctionBegin;
1543bc6112feSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
1544bc6112feSHong Zhang   PetscValidIntPointer(ival,3);
1545bc6112feSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetIcntl_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
1546bc6112feSHong Zhang   PetscFunctionReturn(0);
1547bc6112feSHong Zhang }
1548bc6112feSHong Zhang 
15498928b65cSHong Zhang /* -------------------------------------------------------------------------------------------*/
15508928b65cSHong Zhang #undef __FUNCT__
15518928b65cSHong Zhang #define __FUNCT__ "MatMumpsSetCntl_MUMPS"
15528928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal val)
15538928b65cSHong Zhang {
15548928b65cSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
15558928b65cSHong Zhang 
15568928b65cSHong Zhang   PetscFunctionBegin;
15578928b65cSHong Zhang   mumps->id.CNTL(icntl) = val;
15588928b65cSHong Zhang   PetscFunctionReturn(0);
15598928b65cSHong Zhang }
15608928b65cSHong Zhang 
15618928b65cSHong Zhang #undef __FUNCT__
1562bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetCntl_MUMPS"
1563bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal *val)
1564bc6112feSHong Zhang {
1565bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1566bc6112feSHong Zhang 
1567bc6112feSHong Zhang   PetscFunctionBegin;
1568bc6112feSHong Zhang   *val = mumps->id.CNTL(icntl);
1569bc6112feSHong Zhang   PetscFunctionReturn(0);
1570bc6112feSHong Zhang }
1571bc6112feSHong Zhang 
1572bc6112feSHong Zhang #undef __FUNCT__
15738928b65cSHong Zhang #define __FUNCT__ "MatMumpsSetCntl"
15748928b65cSHong Zhang /*@
15758928b65cSHong Zhang   MatMumpsSetCntl - Set MUMPS parameter CNTL()
15768928b65cSHong Zhang 
15778928b65cSHong Zhang    Logically Collective on Mat
15788928b65cSHong Zhang 
15798928b65cSHong Zhang    Input Parameters:
15808928b65cSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
15818928b65cSHong Zhang .  icntl - index of MUMPS parameter array CNTL()
15828928b65cSHong Zhang -  val - value of MUMPS CNTL(icntl)
15838928b65cSHong Zhang 
15848928b65cSHong Zhang   Options Database:
15858928b65cSHong Zhang .   -mat_mumps_cntl_<icntl> <val>
15868928b65cSHong Zhang 
15878928b65cSHong Zhang    Level: beginner
15888928b65cSHong Zhang 
15898928b65cSHong Zhang    References: MUMPS Users' Guide
15908928b65cSHong Zhang 
15918928b65cSHong Zhang .seealso: MatGetFactor()
15928928b65cSHong Zhang @*/
15938928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl(Mat F,PetscInt icntl,PetscReal val)
15948928b65cSHong Zhang {
15958928b65cSHong Zhang   PetscErrorCode ierr;
15968928b65cSHong Zhang 
15978928b65cSHong Zhang   PetscFunctionBegin;
15988928b65cSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
1599bc6112feSHong Zhang   PetscValidLogicalCollectiveReal(F,val,3);
16008928b65cSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsSetCntl_C",(Mat,PetscInt,PetscReal),(F,icntl,val));CHKERRQ(ierr);
16018928b65cSHong Zhang   PetscFunctionReturn(0);
16028928b65cSHong Zhang }
16038928b65cSHong Zhang 
1604bc6112feSHong Zhang #undef __FUNCT__
1605bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetCntl"
1606a21f80fcSHong Zhang /*@
1607a21f80fcSHong Zhang   MatMumpsGetCntl - Get MUMPS parameter CNTL()
1608a21f80fcSHong Zhang 
1609a21f80fcSHong Zhang    Logically Collective on Mat
1610a21f80fcSHong Zhang 
1611a21f80fcSHong Zhang    Input Parameters:
1612a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1613a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array CNTL()
1614a21f80fcSHong Zhang 
1615a21f80fcSHong Zhang   Output Parameter:
1616a21f80fcSHong Zhang .  val - value of MUMPS CNTL(icntl)
1617a21f80fcSHong Zhang 
1618a21f80fcSHong Zhang    Level: beginner
1619a21f80fcSHong Zhang 
1620a21f80fcSHong Zhang    References: MUMPS Users' Guide
1621a21f80fcSHong Zhang 
1622a21f80fcSHong Zhang .seealso: MatGetFactor()
1623a21f80fcSHong Zhang @*/
1624bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl(Mat F,PetscInt icntl,PetscReal *val)
1625bc6112feSHong Zhang {
1626bc6112feSHong Zhang   PetscErrorCode ierr;
1627bc6112feSHong Zhang 
1628bc6112feSHong Zhang   PetscFunctionBegin;
1629bc6112feSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
1630bc6112feSHong Zhang   PetscValidRealPointer(val,3);
1631bc6112feSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetCntl_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
1632bc6112feSHong Zhang   PetscFunctionReturn(0);
1633bc6112feSHong Zhang }
1634bc6112feSHong Zhang 
1635bc6112feSHong Zhang #undef __FUNCT__
1636ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfo_MUMPS"
1637ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F,PetscInt icntl,PetscInt *info)
1638bc6112feSHong Zhang {
1639bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1640bc6112feSHong Zhang 
1641bc6112feSHong Zhang   PetscFunctionBegin;
1642bc6112feSHong Zhang   *info = mumps->id.INFO(icntl);
1643bc6112feSHong Zhang   PetscFunctionReturn(0);
1644bc6112feSHong Zhang }
1645bc6112feSHong Zhang 
1646bc6112feSHong Zhang #undef __FUNCT__
1647ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfog_MUMPS"
1648ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F,PetscInt icntl,PetscInt *infog)
1649bc6112feSHong Zhang {
1650bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1651bc6112feSHong Zhang 
1652bc6112feSHong Zhang   PetscFunctionBegin;
1653bc6112feSHong Zhang   *infog = mumps->id.INFOG(icntl);
1654bc6112feSHong Zhang   PetscFunctionReturn(0);
1655bc6112feSHong Zhang }
1656bc6112feSHong Zhang 
1657bc6112feSHong Zhang #undef __FUNCT__
1658ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfo_MUMPS"
1659ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfo)
1660bc6112feSHong Zhang {
1661bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1662bc6112feSHong Zhang 
1663bc6112feSHong Zhang   PetscFunctionBegin;
1664bc6112feSHong Zhang   *rinfo = mumps->id.RINFO(icntl);
1665bc6112feSHong Zhang   PetscFunctionReturn(0);
1666bc6112feSHong Zhang }
1667bc6112feSHong Zhang 
1668bc6112feSHong Zhang #undef __FUNCT__
1669ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfog_MUMPS"
1670ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfog)
1671bc6112feSHong Zhang {
1672bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1673bc6112feSHong Zhang 
1674bc6112feSHong Zhang   PetscFunctionBegin;
1675bc6112feSHong Zhang   *rinfog = mumps->id.RINFOG(icntl);
1676bc6112feSHong Zhang   PetscFunctionReturn(0);
1677bc6112feSHong Zhang }
1678bc6112feSHong Zhang 
1679bc6112feSHong Zhang #undef __FUNCT__
1680ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfo"
1681a21f80fcSHong Zhang /*@
1682a21f80fcSHong Zhang   MatMumpsGetInfo - Get MUMPS parameter INFO()
1683a21f80fcSHong Zhang 
1684a21f80fcSHong Zhang    Logically Collective on Mat
1685a21f80fcSHong Zhang 
1686a21f80fcSHong Zhang    Input Parameters:
1687a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1688a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array INFO()
1689a21f80fcSHong Zhang 
1690a21f80fcSHong Zhang   Output Parameter:
1691a21f80fcSHong Zhang .  ival - value of MUMPS INFO(icntl)
1692a21f80fcSHong Zhang 
1693a21f80fcSHong Zhang    Level: beginner
1694a21f80fcSHong Zhang 
1695a21f80fcSHong Zhang    References: MUMPS Users' Guide
1696a21f80fcSHong Zhang 
1697a21f80fcSHong Zhang .seealso: MatGetFactor()
1698a21f80fcSHong Zhang @*/
1699ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo(Mat F,PetscInt icntl,PetscInt *ival)
1700bc6112feSHong Zhang {
1701bc6112feSHong Zhang   PetscErrorCode ierr;
1702bc6112feSHong Zhang 
1703bc6112feSHong Zhang   PetscFunctionBegin;
1704ca810319SHong Zhang   PetscValidIntPointer(ival,3);
1705ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetInfo_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
1706bc6112feSHong Zhang   PetscFunctionReturn(0);
1707bc6112feSHong Zhang }
1708bc6112feSHong Zhang 
1709bc6112feSHong Zhang #undef __FUNCT__
1710ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfog"
1711a21f80fcSHong Zhang /*@
1712a21f80fcSHong Zhang   MatMumpsGetInfog - Get MUMPS parameter INFOG()
1713a21f80fcSHong Zhang 
1714a21f80fcSHong Zhang    Logically Collective on Mat
1715a21f80fcSHong Zhang 
1716a21f80fcSHong Zhang    Input Parameters:
1717a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1718a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array INFOG()
1719a21f80fcSHong Zhang 
1720a21f80fcSHong Zhang   Output Parameter:
1721a21f80fcSHong Zhang .  ival - value of MUMPS INFOG(icntl)
1722a21f80fcSHong Zhang 
1723a21f80fcSHong Zhang    Level: beginner
1724a21f80fcSHong Zhang 
1725a21f80fcSHong Zhang    References: MUMPS Users' Guide
1726a21f80fcSHong Zhang 
1727a21f80fcSHong Zhang .seealso: MatGetFactor()
1728a21f80fcSHong Zhang @*/
1729ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog(Mat F,PetscInt icntl,PetscInt *ival)
1730bc6112feSHong Zhang {
1731bc6112feSHong Zhang   PetscErrorCode ierr;
1732bc6112feSHong Zhang 
1733bc6112feSHong Zhang   PetscFunctionBegin;
1734ca810319SHong Zhang   PetscValidIntPointer(ival,3);
1735ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetInfog_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
1736bc6112feSHong Zhang   PetscFunctionReturn(0);
1737bc6112feSHong Zhang }
1738bc6112feSHong Zhang 
1739bc6112feSHong Zhang #undef __FUNCT__
1740ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfo"
1741a21f80fcSHong Zhang /*@
1742a21f80fcSHong Zhang   MatMumpsGetRinfo - Get MUMPS parameter RINFO()
1743a21f80fcSHong Zhang 
1744a21f80fcSHong Zhang    Logically Collective on Mat
1745a21f80fcSHong Zhang 
1746a21f80fcSHong Zhang    Input Parameters:
1747a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1748a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array RINFO()
1749a21f80fcSHong Zhang 
1750a21f80fcSHong Zhang   Output Parameter:
1751a21f80fcSHong Zhang .  val - value of MUMPS RINFO(icntl)
1752a21f80fcSHong Zhang 
1753a21f80fcSHong Zhang    Level: beginner
1754a21f80fcSHong Zhang 
1755a21f80fcSHong Zhang    References: MUMPS Users' Guide
1756a21f80fcSHong Zhang 
1757a21f80fcSHong Zhang .seealso: MatGetFactor()
1758a21f80fcSHong Zhang @*/
1759ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo(Mat F,PetscInt icntl,PetscReal *val)
1760bc6112feSHong Zhang {
1761bc6112feSHong Zhang   PetscErrorCode ierr;
1762bc6112feSHong Zhang 
1763bc6112feSHong Zhang   PetscFunctionBegin;
1764bc6112feSHong Zhang   PetscValidRealPointer(val,3);
1765ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetRinfo_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
1766bc6112feSHong Zhang   PetscFunctionReturn(0);
1767bc6112feSHong Zhang }
1768bc6112feSHong Zhang 
1769bc6112feSHong Zhang #undef __FUNCT__
1770ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfog"
1771a21f80fcSHong Zhang /*@
1772a21f80fcSHong Zhang   MatMumpsGetRinfog - Get MUMPS parameter RINFOG()
1773a21f80fcSHong Zhang 
1774a21f80fcSHong Zhang    Logically Collective on Mat
1775a21f80fcSHong Zhang 
1776a21f80fcSHong Zhang    Input Parameters:
1777a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1778a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array RINFOG()
1779a21f80fcSHong Zhang 
1780a21f80fcSHong Zhang   Output Parameter:
1781a21f80fcSHong Zhang .  val - value of MUMPS RINFOG(icntl)
1782a21f80fcSHong Zhang 
1783a21f80fcSHong Zhang    Level: beginner
1784a21f80fcSHong Zhang 
1785a21f80fcSHong Zhang    References: MUMPS Users' Guide
1786a21f80fcSHong Zhang 
1787a21f80fcSHong Zhang .seealso: MatGetFactor()
1788a21f80fcSHong Zhang @*/
1789ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog(Mat F,PetscInt icntl,PetscReal *val)
1790bc6112feSHong Zhang {
1791bc6112feSHong Zhang   PetscErrorCode ierr;
1792bc6112feSHong Zhang 
1793bc6112feSHong Zhang   PetscFunctionBegin;
1794bc6112feSHong Zhang   PetscValidRealPointer(val,3);
1795ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetRinfog_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
1796bc6112feSHong Zhang   PetscFunctionReturn(0);
1797bc6112feSHong Zhang }
1798bc6112feSHong Zhang 
179924b6179bSKris Buschelman /*MC
18002692d6eeSBarry Smith   MATSOLVERMUMPS -  A matrix type providing direct solvers (LU and Cholesky) for
180124b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
180224b6179bSKris Buschelman 
180341c8de11SBarry Smith   Works with MATAIJ and MATSBAIJ matrices
180424b6179bSKris Buschelman 
180524b6179bSKris Buschelman   Options Database Keys:
18064e34a73bSHong Zhang +  -mat_mumps_icntl_1 <6>: ICNTL(1): output stream for error messages (None)
18074e34a73bSHong Zhang .  -mat_mumps_icntl_2 <0>: ICNTL(2): output stream for diagnostic printing, statistics, and warning (None)
18084e34a73bSHong Zhang .  -mat_mumps_icntl_3 <0>: ICNTL(3): output stream for global information, collected on the host (None)
18094e34a73bSHong Zhang .  -mat_mumps_icntl_4 <0>: ICNTL(4): level of printing (0 to 4) (None)
18104e34a73bSHong Zhang .  -mat_mumps_icntl_6 <7>: ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7) (None)
18114e34a73bSHong 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)
18124e34a73bSHong Zhang .  -mat_mumps_icntl_8 <77>: ICNTL(8): scaling strategy (-2 to 8 or 77) (None)
18134e34a73bSHong Zhang .  -mat_mumps_icntl_10 <0>: ICNTL(10): max num of refinements (None)
18144e34a73bSHong Zhang .  -mat_mumps_icntl_11 <0>: ICNTL(11): statistics related to an error analysis (via -ksp_view) (None)
18154e34a73bSHong Zhang .  -mat_mumps_icntl_12 <1>: ICNTL(12): an ordering strategy for symmetric matrices (0 to 3) (None)
18164e34a73bSHong Zhang .  -mat_mumps_icntl_13 <0>: ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting (None)
18174e34a73bSHong Zhang .  -mat_mumps_icntl_14 <20>: ICNTL(14): percentage increase in the estimated working space (None)
18184e34a73bSHong Zhang .  -mat_mumps_icntl_19 <0>: ICNTL(19): computes the Schur complement (None)
18194e34a73bSHong Zhang .  -mat_mumps_icntl_22 <0>: ICNTL(22): in-core/out-of-core factorization and solve (0 or 1) (None)
18204e34a73bSHong Zhang .  -mat_mumps_icntl_23 <0>: ICNTL(23): max size of the working memory (MB) that can allocate per processor (None)
18214e34a73bSHong Zhang .  -mat_mumps_icntl_24 <0>: ICNTL(24): detection of null pivot rows (0 or 1) (None)
18224e34a73bSHong Zhang .  -mat_mumps_icntl_25 <0>: ICNTL(25): compute a solution of a deficient matrix and a null space basis (None)
18234e34a73bSHong Zhang .  -mat_mumps_icntl_26 <0>: ICNTL(26): drives the solution phase if a Schur complement matrix (None)
18244e34a73bSHong 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)
18254e34a73bSHong Zhang .  -mat_mumps_icntl_29 <0>: ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis (None)
18264e34a73bSHong Zhang .  -mat_mumps_icntl_30 <0>: ICNTL(30): compute user-specified set of entries in inv(A) (None)
18274e34a73bSHong Zhang .  -mat_mumps_icntl_31 <0>: ICNTL(31): indicates which factors may be discarded during factorization (None)
18284e34a73bSHong Zhang .  -mat_mumps_icntl_33 <0>: ICNTL(33): compute determinant (None)
18294e34a73bSHong Zhang .  -mat_mumps_cntl_1 <0.01>: CNTL(1): relative pivoting threshold (None)
18304e34a73bSHong Zhang .  -mat_mumps_cntl_2 <1.49012e-08>: CNTL(2): stopping criterion of refinement (None)
18314e34a73bSHong Zhang .  -mat_mumps_cntl_3 <0>: CNTL(3): absolute pivoting threshold (None)
18324e34a73bSHong Zhang .  -mat_mumps_cntl_4 <-1>: CNTL(4): value for static pivoting (None)
18334e34a73bSHong Zhang -  -mat_mumps_cntl_5 <0>: CNTL(5): fixation for null pivots (None)
183424b6179bSKris Buschelman 
183524b6179bSKris Buschelman   Level: beginner
183624b6179bSKris Buschelman 
183741c8de11SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage
183841c8de11SBarry Smith 
183924b6179bSKris Buschelman M*/
184024b6179bSKris Buschelman 
184135bd34faSBarry Smith #undef __FUNCT__
184235bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_mumps"
1843f7a08781SBarry Smith static PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type)
184435bd34faSBarry Smith {
184535bd34faSBarry Smith   PetscFunctionBegin;
18462692d6eeSBarry Smith   *type = MATSOLVERMUMPS;
184735bd34faSBarry Smith   PetscFunctionReturn(0);
184835bd34faSBarry Smith }
184935bd34faSBarry Smith 
1850bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI AIJ matrices */
18512877fffaSHong Zhang #undef __FUNCT__
1852bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_aij_mumps"
18538cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
18542877fffaSHong Zhang {
18552877fffaSHong Zhang   Mat            B;
18562877fffaSHong Zhang   PetscErrorCode ierr;
18572877fffaSHong Zhang   Mat_MUMPS      *mumps;
1858ace3abfcSBarry Smith   PetscBool      isSeqAIJ;
18592877fffaSHong Zhang 
18602877fffaSHong Zhang   PetscFunctionBegin;
18612877fffaSHong Zhang   /* Create the factorization matrix */
1862251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr);
1863ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
18642877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
18652877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1866bccb9932SShri Abhyankar   if (isSeqAIJ) {
18670298fd71SBarry Smith     ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr);
1868bccb9932SShri Abhyankar   } else {
18690298fd71SBarry Smith     ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr);
1870bccb9932SShri Abhyankar   }
18712877fffaSHong Zhang 
1872b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
18732205254eSKarl Rupp 
18742877fffaSHong Zhang   B->ops->view        = MatView_MUMPS;
187535bd34faSBarry Smith   B->ops->getinfo     = MatGetInfo_MUMPS;
187620be8e61SHong Zhang   B->ops->getdiagonal = MatGetDiagonal_MUMPS;
18772205254eSKarl Rupp 
1878bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
1879bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
1880bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
1881bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
1882bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
1883bc6112feSHong Zhang 
1884ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
1885ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
1886ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
1887ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
1888450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
1889450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
1890d5f3da31SBarry Smith     B->factortype            = MAT_FACTOR_LU;
1891bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
1892bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
1893746480a1SHong Zhang     mumps->sym = 0;
1894dcd589f8SShri Abhyankar   } else {
189567877ebaSShri Abhyankar     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
1896450b117fSShri Abhyankar     B->factortype                  = MAT_FACTOR_CHOLESKY;
1897bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
1898bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
18996fdc2a6dSBarry Smith     if (A->spd_set && A->spd) mumps->sym = 1;
19006fdc2a6dSBarry Smith     else                      mumps->sym = 2;
1901450b117fSShri Abhyankar   }
19022877fffaSHong Zhang 
19032877fffaSHong Zhang   mumps->isAIJ    = PETSC_TRUE;
1904bf0cc555SLisandro Dalcin   mumps->Destroy  = B->ops->destroy;
19052877fffaSHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
19062877fffaSHong Zhang   B->spptr        = (void*)mumps;
19072205254eSKarl Rupp 
1908f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1909746480a1SHong Zhang 
19102877fffaSHong Zhang   *F = B;
19112877fffaSHong Zhang   PetscFunctionReturn(0);
19122877fffaSHong Zhang }
19132877fffaSHong Zhang 
1914bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI SBAIJ matrices */
19152877fffaSHong Zhang #undef __FUNCT__
1916bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_sbaij_mumps"
19178cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
19182877fffaSHong Zhang {
19192877fffaSHong Zhang   Mat            B;
19202877fffaSHong Zhang   PetscErrorCode ierr;
19212877fffaSHong Zhang   Mat_MUMPS      *mumps;
1922ace3abfcSBarry Smith   PetscBool      isSeqSBAIJ;
19232877fffaSHong Zhang 
19242877fffaSHong Zhang   PetscFunctionBegin;
1925ce94432eSBarry Smith   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
1926ce94432eSBarry 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");
1927251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
19282877fffaSHong Zhang   /* Create the factorization matrix */
1929ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
19302877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
19312877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1932b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
1933bccb9932SShri Abhyankar   if (isSeqSBAIJ) {
19340298fd71SBarry Smith     ierr = MatSeqSBAIJSetPreallocation(B,1,0,NULL);CHKERRQ(ierr);
19352205254eSKarl Rupp 
193616ebf90aSShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
1937dcd589f8SShri Abhyankar   } else {
19380298fd71SBarry Smith     ierr = MatMPISBAIJSetPreallocation(B,1,0,NULL,0,NULL);CHKERRQ(ierr);
19392205254eSKarl Rupp 
1940bccb9932SShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
1941bccb9932SShri Abhyankar   }
1942bccb9932SShri Abhyankar 
194367877ebaSShri Abhyankar   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
1944bccb9932SShri Abhyankar   B->ops->view                   = MatView_MUMPS;
194520be8e61SHong Zhang   B->ops->getdiagonal            = MatGetDiagonal_MUMPS;
19462205254eSKarl Rupp 
1947bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
1948b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
1949b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
1950b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
1951b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
1952bc6112feSHong Zhang 
1953ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
1954ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
1955ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
1956ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
19572205254eSKarl Rupp 
1958f4762488SHong Zhang   B->factortype = MAT_FACTOR_CHOLESKY;
19596fdc2a6dSBarry Smith   if (A->spd_set && A->spd) mumps->sym = 1;
19606fdc2a6dSBarry Smith   else                      mumps->sym = 2;
1961a214ac2aSShri Abhyankar 
1962bccb9932SShri Abhyankar   mumps->isAIJ    = PETSC_FALSE;
1963bf0cc555SLisandro Dalcin   mumps->Destroy  = B->ops->destroy;
1964f3c0ef26SHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
19652877fffaSHong Zhang   B->spptr        = (void*)mumps;
19662205254eSKarl Rupp 
1967f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1968746480a1SHong Zhang 
19692877fffaSHong Zhang   *F = B;
19702877fffaSHong Zhang   PetscFunctionReturn(0);
19712877fffaSHong Zhang }
197297969023SHong Zhang 
1973450b117fSShri Abhyankar #undef __FUNCT__
1974bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_baij_mumps"
19758cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
197667877ebaSShri Abhyankar {
197767877ebaSShri Abhyankar   Mat            B;
197867877ebaSShri Abhyankar   PetscErrorCode ierr;
197967877ebaSShri Abhyankar   Mat_MUMPS      *mumps;
1980ace3abfcSBarry Smith   PetscBool      isSeqBAIJ;
198167877ebaSShri Abhyankar 
198267877ebaSShri Abhyankar   PetscFunctionBegin;
198367877ebaSShri Abhyankar   /* Create the factorization matrix */
1984251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr);
1985ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
198667877ebaSShri Abhyankar   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
198767877ebaSShri Abhyankar   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1988bccb9932SShri Abhyankar   if (isSeqBAIJ) {
19890298fd71SBarry Smith     ierr = MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,NULL);CHKERRQ(ierr);
1990bccb9932SShri Abhyankar   } else {
19910298fd71SBarry Smith     ierr = MatMPIBAIJSetPreallocation(B,A->rmap->bs,0,NULL,0,NULL);CHKERRQ(ierr);
1992bccb9932SShri Abhyankar   }
1993450b117fSShri Abhyankar 
1994b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
1995450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
1996450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
1997450b117fSShri Abhyankar     B->factortype            = MAT_FACTOR_LU;
1998bccb9932SShri Abhyankar     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
1999bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
2000746480a1SHong Zhang     mumps->sym = 0;
2001f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n");
2002bccb9932SShri Abhyankar 
2003450b117fSShri Abhyankar   B->ops->view        = MatView_MUMPS;
200420be8e61SHong Zhang   B->ops->getdiagonal = MatGetDiagonal_MUMPS;
20052205254eSKarl Rupp 
2006bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
2007bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
2008bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
2009bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
2010bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
2011bc6112feSHong Zhang 
2012ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
2013ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
2014ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
2015ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
2016450b117fSShri Abhyankar 
2017450b117fSShri Abhyankar   mumps->isAIJ    = PETSC_TRUE;
2018bf0cc555SLisandro Dalcin   mumps->Destroy  = B->ops->destroy;
2019450b117fSShri Abhyankar   B->ops->destroy = MatDestroy_MUMPS;
2020450b117fSShri Abhyankar   B->spptr        = (void*)mumps;
20212205254eSKarl Rupp 
2022f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
2023746480a1SHong Zhang 
2024450b117fSShri Abhyankar   *F = B;
2025450b117fSShri Abhyankar   PetscFunctionReturn(0);
2026450b117fSShri Abhyankar }
202742c9c57cSBarry Smith 
202842c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mumps(Mat,MatFactorType,Mat*);
202942c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_baij_mumps(Mat,MatFactorType,Mat*);
203042c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat,MatFactorType,Mat*);
203142c9c57cSBarry Smith 
203242c9c57cSBarry Smith #undef __FUNCT__
203342c9c57cSBarry Smith #define __FUNCT__ "MatSolverPackageRegister_MUMPS"
203429b38603SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_MUMPS(void)
203542c9c57cSBarry Smith {
203642c9c57cSBarry Smith   PetscErrorCode ierr;
203742c9c57cSBarry Smith 
203842c9c57cSBarry Smith   PetscFunctionBegin;
203942c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ,         MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr);
204042c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ,         MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr);
204142c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ,         MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr);
204242c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ,         MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr);
204342c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPISBAIJ,         MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr);
204442c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ,         MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr);
204542c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ,         MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr);
204642c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ,         MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr);
204742c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ,         MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr);
204842c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQSBAIJ,         MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr);
204942c9c57cSBarry Smith   PetscFunctionReturn(0);
205042c9c57cSBarry Smith }
205142c9c57cSBarry Smith 
2052