xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision d341cd04ead25786cb579bd428dce3b9509eec8a)
11c2a3de1SBarry Smith 
2397b6df1SKris Buschelman /*
3c2b5dc30SHong Zhang     Provides an interface to the MUMPS sparse solver
4397b6df1SKris Buschelman */
551d5961aSHong Zhang 
6c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> /*I  "petscmat.h"  I*/
7c6db04a5SJed Brown #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
8397b6df1SKris Buschelman 
9397b6df1SKris Buschelman EXTERN_C_BEGIN
10397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
112907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
122907cef9SHong Zhang #include <cmumps_c.h>
132907cef9SHong Zhang #else
14c6db04a5SJed Brown #include <zmumps_c.h>
152907cef9SHong Zhang #endif
162907cef9SHong Zhang #else
172907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
182907cef9SHong Zhang #include <smumps_c.h>
19397b6df1SKris Buschelman #else
20c6db04a5SJed Brown #include <dmumps_c.h>
21397b6df1SKris Buschelman #endif
222907cef9SHong Zhang #endif
23397b6df1SKris Buschelman EXTERN_C_END
24397b6df1SKris Buschelman #define JOB_INIT -1
253d472b54SHong Zhang #define JOB_FACTSYMBOLIC 1
263d472b54SHong Zhang #define JOB_FACTNUMERIC 2
273d472b54SHong Zhang #define JOB_SOLVE 3
28397b6df1SKris Buschelman #define JOB_END -2
293d472b54SHong Zhang 
302907cef9SHong Zhang /* calls to MUMPS */
312907cef9SHong Zhang #if defined(PETSC_USE_COMPLEX)
322907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
332907cef9SHong Zhang #define PetscMUMPS_c cmumps_c
342907cef9SHong Zhang #else
352907cef9SHong Zhang #define PetscMUMPS_c zmumps_c
362907cef9SHong Zhang #endif
372907cef9SHong Zhang #else
382907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
392907cef9SHong Zhang #define PetscMUMPS_c smumps_c
402907cef9SHong Zhang #else
412907cef9SHong Zhang #define PetscMUMPS_c dmumps_c
422907cef9SHong Zhang #endif
432907cef9SHong Zhang #endif
442907cef9SHong Zhang 
453d472b54SHong Zhang 
46397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */
47397b6df1SKris Buschelman #define ICNTL(I) icntl[(I)-1]
48397b6df1SKris Buschelman #define CNTL(I) cntl[(I)-1]
49397b6df1SKris Buschelman #define INFOG(I) infog[(I)-1]
50a7aca84bSHong Zhang #define INFO(I) info[(I)-1]
51397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1]
52adc1d99fSHong Zhang #define RINFO(I) rinfo[(I)-1]
53397b6df1SKris Buschelman 
54397b6df1SKris Buschelman typedef struct {
55397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
562907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
572907cef9SHong Zhang   CMUMPS_STRUC_C id;
582907cef9SHong Zhang #else
59397b6df1SKris Buschelman   ZMUMPS_STRUC_C id;
602907cef9SHong Zhang #endif
612907cef9SHong Zhang #else
622907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
632907cef9SHong Zhang   SMUMPS_STRUC_C id;
64397b6df1SKris Buschelman #else
65397b6df1SKris Buschelman   DMUMPS_STRUC_C id;
66397b6df1SKris Buschelman #endif
672907cef9SHong Zhang #endif
682907cef9SHong Zhang 
69397b6df1SKris Buschelman   MatStructure matstruc;
70c1490034SHong Zhang   PetscMPIInt  myid,size;
71a5e57a09SHong Zhang   PetscInt     *irn,*jcn,nz,sym;
72397b6df1SKris Buschelman   PetscScalar  *val;
73397b6df1SKris Buschelman   MPI_Comm     comm_mumps;
74329ec9b3SHong Zhang   VecScatter   scat_rhs, scat_sol;
7564e6c443SBarry Smith   PetscBool    isAIJ,CleanUpMUMPS;
76329ec9b3SHong Zhang   Vec          b_seq,x_seq;
77a5e57a09SHong Zhang   PetscInt     ICNTL9_pre;   /* check if ICNTL(9) is changed from previous MatSolve */
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 /*
86*d341cd04SHong Zhang   MatConvertToTriples_A_B - convert Petsc matrix to triples: row[nz], col[nz], val[nz]
87*d341cd04SHong 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 = VecDestroy(&mumps->b_seq);CHKERRQ(ierr);
547a5e57a09SHong Zhang     ierr = VecScatterDestroy(&mumps->scat_sol);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;
6532205254eSKarl Rupp 
6540ad0caddSJed Brown   ierr = MatSolve_MUMPS(A,b,x);CHKERRQ(ierr);
6552205254eSKarl Rupp 
656a5e57a09SHong Zhang   mumps->id.ICNTL(9) = 1;
65751d5961aSHong Zhang   PetscFunctionReturn(0);
65851d5961aSHong Zhang }
65951d5961aSHong Zhang 
660e0b74bf9SHong Zhang #undef __FUNCT__
661e0b74bf9SHong Zhang #define __FUNCT__ "MatMatSolve_MUMPS"
662e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X)
663e0b74bf9SHong Zhang {
664bda8bf91SBarry Smith   PetscErrorCode ierr;
665bda8bf91SBarry Smith   PetscBool      flg;
666bda8bf91SBarry Smith 
667e0b74bf9SHong Zhang   PetscFunctionBegin;
6680298fd71SBarry Smith   ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
669ce94432eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
6700298fd71SBarry Smith   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
671ce94432eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
6722205254eSKarl Rupp   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_MUMPS() is not implemented yet");
673e0b74bf9SHong Zhang   PetscFunctionReturn(0);
674e0b74bf9SHong Zhang }
675e0b74bf9SHong Zhang 
676ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX)
677a58c3f20SHong Zhang /*
678a58c3f20SHong Zhang   input:
679a58c3f20SHong Zhang    F:        numeric factor
680a58c3f20SHong Zhang   output:
681a58c3f20SHong Zhang    nneg:     total number of negative pivots
682a58c3f20SHong Zhang    nzero:    0
683a58c3f20SHong Zhang    npos:     (global dimension of F) - nneg
684a58c3f20SHong Zhang */
685a58c3f20SHong Zhang 
686a58c3f20SHong Zhang #undef __FUNCT__
687a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
688dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
689a58c3f20SHong Zhang {
690a5e57a09SHong Zhang   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
691dfbe8321SBarry Smith   PetscErrorCode ierr;
692c1490034SHong Zhang   PetscMPIInt    size;
693a58c3f20SHong Zhang 
694a58c3f20SHong Zhang   PetscFunctionBegin;
695ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)F),&size);CHKERRQ(ierr);
696bcb30aebSHong 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 */
697a5e57a09SHong 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));
698ed85ac9fSHong Zhang 
699710ac8efSHong Zhang   if (nneg) *nneg = mumps->id.INFOG(12);
700ed85ac9fSHong Zhang   if (nzero || npos) {
701ed85ac9fSHong 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");
702710ac8efSHong Zhang     if (nzero) *nzero = mumps->id.INFOG(28);
703710ac8efSHong Zhang     if (npos) *npos   = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28));
704a58c3f20SHong Zhang   }
705a58c3f20SHong Zhang   PetscFunctionReturn(0);
706a58c3f20SHong Zhang }
707ace3df97SHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */
708a58c3f20SHong Zhang 
709397b6df1SKris Buschelman #undef __FUNCT__
710f6c57405SHong Zhang #define __FUNCT__ "MatFactorNumeric_MUMPS"
7110481f469SBarry Smith PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
712af281ebdSHong Zhang {
713a5e57a09SHong Zhang   Mat_MUMPS      *mumps =(Mat_MUMPS*)(F)->spptr;
7146849ba73SBarry Smith   PetscErrorCode ierr;
715e09efc27SHong Zhang   Mat            F_diag;
716ace3abfcSBarry Smith   PetscBool      isMPIAIJ;
717397b6df1SKris Buschelman 
718397b6df1SKris Buschelman   PetscFunctionBegin;
719a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
720397b6df1SKris Buschelman 
721397b6df1SKris Buschelman   /* numerical factorization phase */
722329ec9b3SHong Zhang   /*-------------------------------*/
723a5e57a09SHong Zhang   mumps->id.job = JOB_FACTNUMERIC;
724a5e57a09SHong Zhang   if (!mumps->id.ICNTL(18)) {
725a5e57a09SHong Zhang     if (!mumps->myid) {
726397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
7272907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
728a5e57a09SHong Zhang       mumps->id.a = (mumps_complex*)mumps->val;
7292907cef9SHong Zhang #else
730a5e57a09SHong Zhang       mumps->id.a = (mumps_double_complex*)mumps->val;
7312907cef9SHong Zhang #endif
732397b6df1SKris Buschelman #else
733a5e57a09SHong Zhang       mumps->id.a = mumps->val;
734397b6df1SKris Buschelman #endif
735397b6df1SKris Buschelman     }
736397b6df1SKris Buschelman   } else {
737397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
7382907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
739a5e57a09SHong Zhang     mumps->id.a_loc = (mumps_complex*)mumps->val;
7402907cef9SHong Zhang #else
741a5e57a09SHong Zhang     mumps->id.a_loc = (mumps_double_complex*)mumps->val;
7422907cef9SHong Zhang #endif
743397b6df1SKris Buschelman #else
744a5e57a09SHong Zhang     mumps->id.a_loc = mumps->val;
745397b6df1SKris Buschelman #endif
746397b6df1SKris Buschelman   }
747a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
748a5e57a09SHong Zhang   if (mumps->id.INFOG(1) < 0) {
749151787a6SHong Zhang     if (mumps->id.INFO(1) == -13) {
750151787a6SHong Zhang       if (mumps->id.INFO(2) < 0) {
751151787a6SHong 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));
752151787a6SHong Zhang       } else {
753151787a6SHong 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));
754151787a6SHong Zhang       }
755151787a6SHong 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));
756397b6df1SKris Buschelman   }
757a5e57a09SHong 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));
758397b6df1SKris Buschelman 
759dcd589f8SShri Abhyankar   (F)->assembled      = PETSC_TRUE;
760a5e57a09SHong Zhang   mumps->matstruc     = SAME_NONZERO_PATTERN;
761a5e57a09SHong Zhang   mumps->CleanUpMUMPS = PETSC_TRUE;
76267877ebaSShri Abhyankar 
763a5e57a09SHong Zhang   if (mumps->size > 1) {
76467877ebaSShri Abhyankar     PetscInt    lsol_loc;
76567877ebaSShri Abhyankar     PetscScalar *sol_loc;
7662205254eSKarl Rupp 
767c2093ab7SHong Zhang     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr);
768c2093ab7SHong Zhang     if (isMPIAIJ) F_diag = ((Mat_MPIAIJ*)(F)->data)->A;
769c2093ab7SHong Zhang     else F_diag = ((Mat_MPISBAIJ*)(F)->data)->A;
770c2093ab7SHong Zhang     F_diag->assembled = PETSC_TRUE;
771c2093ab7SHong Zhang 
772c2093ab7SHong Zhang     /* distributed solution; Create x_seq=sol_loc for repeated use */
773c2093ab7SHong Zhang     if (mumps->x_seq) {
774c2093ab7SHong Zhang       ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
775c2093ab7SHong Zhang       ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr);
776c2093ab7SHong Zhang       ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr);
777c2093ab7SHong Zhang     }
778a5e57a09SHong Zhang     lsol_loc = mumps->id.INFO(23); /* length of sol_loc */
779dcca6d9dSJed Brown     ierr = PetscMalloc2(lsol_loc,&sol_loc,lsol_loc,&mumps->id.isol_loc);CHKERRQ(ierr);
780a5e57a09SHong Zhang     mumps->id.lsol_loc = lsol_loc;
78167877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
7822907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
783a5e57a09SHong Zhang     mumps->id.sol_loc = (mumps_complex*)sol_loc;
7842907cef9SHong Zhang #else
785a5e57a09SHong Zhang     mumps->id.sol_loc = (mumps_double_complex*)sol_loc;
7862907cef9SHong Zhang #endif
78767877ebaSShri Abhyankar #else
788a5e57a09SHong Zhang     mumps->id.sol_loc = sol_loc;
78967877ebaSShri Abhyankar #endif
790a5e57a09SHong Zhang     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq);CHKERRQ(ierr);
79167877ebaSShri Abhyankar   }
792397b6df1SKris Buschelman   PetscFunctionReturn(0);
793397b6df1SKris Buschelman }
794397b6df1SKris Buschelman 
7959a2535b5SHong Zhang /* Sets MUMPS options from the options database */
796dcd589f8SShri Abhyankar #undef __FUNCT__
7979a2535b5SHong Zhang #define __FUNCT__ "PetscSetMUMPSFromOptions"
7989a2535b5SHong Zhang PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A)
799dcd589f8SShri Abhyankar {
8009a2535b5SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
801dcd589f8SShri Abhyankar   PetscErrorCode ierr;
802dcd589f8SShri Abhyankar   PetscInt       icntl;
803ace3abfcSBarry Smith   PetscBool      flg;
804dcd589f8SShri Abhyankar 
805dcd589f8SShri Abhyankar   PetscFunctionBegin;
806ce94432eSBarry Smith   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
8079a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);CHKERRQ(ierr);
8089a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(1) = icntl;
8099a2535b5SHong 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);
8109a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(2) = icntl;
8119a2535b5SHong 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);
8129a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(3) = icntl;
813dcd589f8SShri Abhyankar 
8149a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",mumps->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
8159a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(4) = icntl;
8169a2535b5SHong Zhang   if (mumps->id.ICNTL(4) || PetscLogPrintInfo) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */
8179a2535b5SHong Zhang 
818*d341cd04SHong 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);
8199a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(6) = icntl;
8209a2535b5SHong Zhang 
821*d341cd04SHong 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);
822dcd589f8SShri Abhyankar   if (flg) {
8232205254eSKarl 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");
8242205254eSKarl Rupp     else mumps->id.ICNTL(7) = icntl;
825dcd589f8SShri Abhyankar   }
826e0b74bf9SHong Zhang 
8270298fd71SBarry 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);
828*d341cd04SHong 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() */
8290298fd71SBarry 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);
830*d341cd04SHong 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);
831*d341cd04SHong 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);
832*d341cd04SHong 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);
833*d341cd04SHong 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);
834*d341cd04SHong 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);
835*d341cd04SHong 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); -- non-dense rhs is not supported in this interface */
836*d341cd04SHong 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 */
8379a2535b5SHong Zhang 
838*d341cd04SHong 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);
8390298fd71SBarry 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);
8400298fd71SBarry 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);
8419a2535b5SHong Zhang   if (mumps->id.ICNTL(24)) {
8429a2535b5SHong Zhang     mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */
843d7ebd59bSHong Zhang   }
844d7ebd59bSHong Zhang 
845*d341cd04SHong 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);
846*d341cd04SHong 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);
847*d341cd04SHong 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); //work on MatSolve_MUMPS()!!!
8480298fd71SBarry 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);
849*d341cd04SHong 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);
8500298fd71SBarry 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);
851*d341cd04SHong 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);
852*d341cd04SHong 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);
8530298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),NULL);CHKERRQ(ierr);
854dcd589f8SShri Abhyankar 
8550298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),NULL);CHKERRQ(ierr);
8560298fd71SBarry 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);
8570298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),NULL);CHKERRQ(ierr);
8580298fd71SBarry 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);
8590298fd71SBarry 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);
860e5bb22a1SHong Zhang 
8610298fd71SBarry Smith   ierr = PetscOptionsString("-mat_mumps_ooc_tmpdir", "out of core directory", "None", mumps->id.ooc_tmpdir, mumps->id.ooc_tmpdir, 256, NULL);
862dcd589f8SShri Abhyankar   PetscOptionsEnd();
863dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
864dcd589f8SShri Abhyankar }
865dcd589f8SShri Abhyankar 
866dcd589f8SShri Abhyankar #undef __FUNCT__
867dcd589f8SShri Abhyankar #define __FUNCT__ "PetscInitializeMUMPS"
868f697e70eSHong Zhang PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps)
869dcd589f8SShri Abhyankar {
870dcd589f8SShri Abhyankar   PetscErrorCode ierr;
871dcd589f8SShri Abhyankar 
872dcd589f8SShri Abhyankar   PetscFunctionBegin;
873ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)A), &mumps->myid);
874ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&mumps->size);CHKERRQ(ierr);
875ce94432eSBarry Smith   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(mumps->comm_mumps));CHKERRQ(ierr);
8762205254eSKarl Rupp 
877f697e70eSHong Zhang   mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps);
878f697e70eSHong Zhang 
879f697e70eSHong Zhang   mumps->id.job = JOB_INIT;
880f697e70eSHong Zhang   mumps->id.par = 1;  /* host participates factorizaton and solve */
881f697e70eSHong Zhang   mumps->id.sym = mumps->sym;
8822907cef9SHong Zhang   PetscMUMPS_c(&mumps->id);
883f697e70eSHong Zhang 
884f697e70eSHong Zhang   mumps->CleanUpMUMPS = PETSC_FALSE;
8850298fd71SBarry Smith   mumps->scat_rhs     = NULL;
8860298fd71SBarry Smith   mumps->scat_sol     = NULL;
8879a2535b5SHong Zhang 
88870544d5fSHong Zhang   /* set PETSc-MUMPS default options - override MUMPS default */
8899a2535b5SHong Zhang   mumps->id.ICNTL(3) = 0;
8909a2535b5SHong Zhang   mumps->id.ICNTL(4) = 0;
8919a2535b5SHong Zhang   if (mumps->size == 1) {
8929a2535b5SHong Zhang     mumps->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
8939a2535b5SHong Zhang   } else {
8949a2535b5SHong Zhang     mumps->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
89570544d5fSHong Zhang     mumps->id.ICNTL(21) = 1;   /* distributed solution */
8969a2535b5SHong Zhang   }
897dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
898dcd589f8SShri Abhyankar }
899dcd589f8SShri Abhyankar 
900a5e57a09SHong 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 */
901397b6df1SKris Buschelman #undef __FUNCT__
902f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
9030481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
904b24902e0SBarry Smith {
905a5e57a09SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
906dcd589f8SShri Abhyankar   PetscErrorCode ierr;
90767877ebaSShri Abhyankar   Vec            b;
90867877ebaSShri Abhyankar   IS             is_iden;
90967877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
910397b6df1SKris Buschelman 
911397b6df1SKris Buschelman   PetscFunctionBegin;
912a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
913dcd589f8SShri Abhyankar 
9149a2535b5SHong Zhang   /* Set MUMPS options from the options database */
9159a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
916dcd589f8SShri Abhyankar 
917a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
918dcd589f8SShri Abhyankar 
91967877ebaSShri Abhyankar   /* analysis phase */
92067877ebaSShri Abhyankar   /*----------------*/
921a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
922a5e57a09SHong Zhang   mumps->id.n   = M;
923a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
92467877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
925a5e57a09SHong Zhang     if (!mumps->myid) {
926a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
927a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
92867877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
9292907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
930a5e57a09SHong Zhang         mumps->id.a = (mumps_complex*)mumps->val;
9312907cef9SHong Zhang #else
932a5e57a09SHong Zhang         mumps->id.a = (mumps_double_complex*)mumps->val;
9332907cef9SHong Zhang #endif
93467877ebaSShri Abhyankar #else
935a5e57a09SHong Zhang         mumps->id.a = mumps->val;
93667877ebaSShri Abhyankar #endif
93767877ebaSShri Abhyankar       }
938a5e57a09SHong Zhang       if (mumps->id.ICNTL(7) == 1) { /* use user-provide matrix ordering - assuming r = c ordering */
9395248a706SHong Zhang         /*
9405248a706SHong Zhang         PetscBool      flag;
9415248a706SHong Zhang         ierr = ISEqual(r,c,&flag);CHKERRQ(ierr);
9425248a706SHong Zhang         if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"row_perm != col_perm");
9435248a706SHong Zhang         ierr = ISView(r,PETSC_VIEWER_STDOUT_SELF);
9445248a706SHong Zhang          */
945a5e57a09SHong Zhang         if (!mumps->myid) {
946e0b74bf9SHong Zhang           const PetscInt *idx;
947e0b74bf9SHong Zhang           PetscInt       i,*perm_in;
9482205254eSKarl Rupp 
949785e854fSJed Brown           ierr = PetscMalloc1(M,&perm_in);CHKERRQ(ierr);
950e0b74bf9SHong Zhang           ierr = ISGetIndices(r,&idx);CHKERRQ(ierr);
9512205254eSKarl Rupp 
952a5e57a09SHong Zhang           mumps->id.perm_in = perm_in;
953e0b74bf9SHong Zhang           for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */
954e0b74bf9SHong Zhang           ierr = ISRestoreIndices(r,&idx);CHKERRQ(ierr);
955e0b74bf9SHong Zhang         }
956e0b74bf9SHong Zhang       }
95767877ebaSShri Abhyankar     }
95867877ebaSShri Abhyankar     break;
95967877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
960a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
961a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
962a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
96367877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
9642907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
965a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_complex*)mumps->val;
9662907cef9SHong Zhang #else
967a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_double_complex*)mumps->val;
9682907cef9SHong Zhang #endif
96967877ebaSShri Abhyankar #else
970a5e57a09SHong Zhang       mumps->id.a_loc = mumps->val;
97167877ebaSShri Abhyankar #endif
97267877ebaSShri Abhyankar     }
97367877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
974a5e57a09SHong Zhang     if (!mumps->myid) {
975a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
97667877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
97767877ebaSShri Abhyankar     } else {
978a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
97967877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
98067877ebaSShri Abhyankar     }
9812a7a6963SBarry Smith     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
982a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
9836bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
9846bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
98567877ebaSShri Abhyankar     break;
98667877ebaSShri Abhyankar   }
987a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
988a5e57a09SHong 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));
98967877ebaSShri Abhyankar 
990719d5645SBarry Smith   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
991dcd589f8SShri Abhyankar   F->ops->solve           = MatSolve_MUMPS;
99251d5961aSHong Zhang   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
99317f96c7aSHong Zhang   F->ops->matsolve        = 0;  /* use MatMatSolve_Basic() until mumps supports distributed rhs */
994b24902e0SBarry Smith   PetscFunctionReturn(0);
995b24902e0SBarry Smith }
996b24902e0SBarry Smith 
997450b117fSShri Abhyankar /* Note the Petsc r and c permutations are ignored */
998450b117fSShri Abhyankar #undef __FUNCT__
999450b117fSShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_BAIJMUMPS"
1000450b117fSShri Abhyankar PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1001450b117fSShri Abhyankar {
1002a5e57a09SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
1003dcd589f8SShri Abhyankar   PetscErrorCode ierr;
100467877ebaSShri Abhyankar   Vec            b;
100567877ebaSShri Abhyankar   IS             is_iden;
100667877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1007450b117fSShri Abhyankar 
1008450b117fSShri Abhyankar   PetscFunctionBegin;
1009a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1010dcd589f8SShri Abhyankar 
10119a2535b5SHong Zhang   /* Set MUMPS options from the options database */
10129a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1013dcd589f8SShri Abhyankar 
1014a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
101567877ebaSShri Abhyankar 
101667877ebaSShri Abhyankar   /* analysis phase */
101767877ebaSShri Abhyankar   /*----------------*/
1018a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1019a5e57a09SHong Zhang   mumps->id.n   = M;
1020a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
102167877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1022a5e57a09SHong Zhang     if (!mumps->myid) {
1023a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1024a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
102567877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
10262907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1027a5e57a09SHong Zhang         mumps->id.a = (mumps_complex*)mumps->val;
10282907cef9SHong Zhang #else
1029a5e57a09SHong Zhang         mumps->id.a = (mumps_double_complex*)mumps->val;
10302907cef9SHong Zhang #endif
103167877ebaSShri Abhyankar #else
1032a5e57a09SHong Zhang         mumps->id.a = mumps->val;
103367877ebaSShri Abhyankar #endif
103467877ebaSShri Abhyankar       }
103567877ebaSShri Abhyankar     }
103667877ebaSShri Abhyankar     break;
103767877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1038a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1039a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1040a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
104167877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
10422907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1043a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_complex*)mumps->val;
10442907cef9SHong Zhang #else
1045a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_double_complex*)mumps->val;
10462907cef9SHong Zhang #endif
104767877ebaSShri Abhyankar #else
1048a5e57a09SHong Zhang       mumps->id.a_loc = mumps->val;
104967877ebaSShri Abhyankar #endif
105067877ebaSShri Abhyankar     }
105167877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1052a5e57a09SHong Zhang     if (!mumps->myid) {
1053a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
105467877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
105567877ebaSShri Abhyankar     } else {
1056a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
105767877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
105867877ebaSShri Abhyankar     }
10592a7a6963SBarry Smith     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1060a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
10616bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
10626bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
106367877ebaSShri Abhyankar     break;
106467877ebaSShri Abhyankar   }
1065a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
1066a5e57a09SHong 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));
106767877ebaSShri Abhyankar 
1068450b117fSShri Abhyankar   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1069dcd589f8SShri Abhyankar   F->ops->solve           = MatSolve_MUMPS;
107051d5961aSHong Zhang   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1071450b117fSShri Abhyankar   PetscFunctionReturn(0);
1072450b117fSShri Abhyankar }
1073b24902e0SBarry Smith 
1074141f4205SHong Zhang /* Note the Petsc r permutation and factor info are ignored */
1075397b6df1SKris Buschelman #undef __FUNCT__
107667877ebaSShri Abhyankar #define __FUNCT__ "MatCholeskyFactorSymbolic_MUMPS"
107767877ebaSShri Abhyankar PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
1078b24902e0SBarry Smith {
1079a5e57a09SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
1080dcd589f8SShri Abhyankar   PetscErrorCode ierr;
108167877ebaSShri Abhyankar   Vec            b;
108267877ebaSShri Abhyankar   IS             is_iden;
108367877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1084397b6df1SKris Buschelman 
1085397b6df1SKris Buschelman   PetscFunctionBegin;
1086a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1087dcd589f8SShri Abhyankar 
10889a2535b5SHong Zhang   /* Set MUMPS options from the options database */
10899a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1090dcd589f8SShri Abhyankar 
1091a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1092dcd589f8SShri Abhyankar 
109367877ebaSShri Abhyankar   /* analysis phase */
109467877ebaSShri Abhyankar   /*----------------*/
1095a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1096a5e57a09SHong Zhang   mumps->id.n   = M;
1097a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
109867877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1099a5e57a09SHong Zhang     if (!mumps->myid) {
1100a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1101a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
110267877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
11032907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1104a5e57a09SHong Zhang         mumps->id.a = (mumps_complex*)mumps->val;
11052907cef9SHong Zhang #else
1106a5e57a09SHong Zhang         mumps->id.a = (mumps_double_complex*)mumps->val;
11072907cef9SHong Zhang #endif
110867877ebaSShri Abhyankar #else
1109a5e57a09SHong Zhang         mumps->id.a = mumps->val;
111067877ebaSShri Abhyankar #endif
111167877ebaSShri Abhyankar       }
111267877ebaSShri Abhyankar     }
111367877ebaSShri Abhyankar     break;
111467877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1115a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1116a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1117a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
111867877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
11192907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1120a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_complex*)mumps->val;
11212907cef9SHong Zhang #else
1122a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_double_complex*)mumps->val;
11232907cef9SHong Zhang #endif
112467877ebaSShri Abhyankar #else
1125a5e57a09SHong Zhang       mumps->id.a_loc = mumps->val;
112667877ebaSShri Abhyankar #endif
112767877ebaSShri Abhyankar     }
112867877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1129a5e57a09SHong Zhang     if (!mumps->myid) {
1130a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
113167877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
113267877ebaSShri Abhyankar     } else {
1133a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
113467877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
113567877ebaSShri Abhyankar     }
11362a7a6963SBarry Smith     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1137a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
11386bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
11396bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
114067877ebaSShri Abhyankar     break;
114167877ebaSShri Abhyankar   }
1142a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
1143a5e57a09SHong 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));
114467877ebaSShri Abhyankar 
11452792810eSHong Zhang   F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
1146dcd589f8SShri Abhyankar   F->ops->solve                 = MatSolve_MUMPS;
114751d5961aSHong Zhang   F->ops->solvetranspose        = MatSolve_MUMPS;
114830c107b7SHong Zhang   F->ops->matsolve              = 0; /* use MatMatSolve_Basic() until mumps supports distributed rhs */
1149db4efbfdSBarry Smith #if !defined(PETSC_USE_COMPLEX)
115005aa0992SJose Roman   F->ops->getinertia = MatGetInertia_SBAIJMUMPS;
115105aa0992SJose Roman #else
11520298fd71SBarry Smith   F->ops->getinertia = NULL;
1153db4efbfdSBarry Smith #endif
1154b24902e0SBarry Smith   PetscFunctionReturn(0);
1155b24902e0SBarry Smith }
1156b24902e0SBarry Smith 
1157397b6df1SKris Buschelman #undef __FUNCT__
115864e6c443SBarry Smith #define __FUNCT__ "MatView_MUMPS"
115964e6c443SBarry Smith PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
116074ed9c26SBarry Smith {
1161f6c57405SHong Zhang   PetscErrorCode    ierr;
116264e6c443SBarry Smith   PetscBool         iascii;
116364e6c443SBarry Smith   PetscViewerFormat format;
1164a5e57a09SHong Zhang   Mat_MUMPS         *mumps=(Mat_MUMPS*)A->spptr;
1165f6c57405SHong Zhang 
1166f6c57405SHong Zhang   PetscFunctionBegin;
116764e6c443SBarry Smith   /* check if matrix is mumps type */
116864e6c443SBarry Smith   if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0);
116964e6c443SBarry Smith 
1170251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
117164e6c443SBarry Smith   if (iascii) {
117264e6c443SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
117364e6c443SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO) {
117464e6c443SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
1175a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                   %d \n",mumps->id.sym);CHKERRQ(ierr);
1176a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):            %d \n",mumps->id.par);CHKERRQ(ierr);
1177a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):         %d \n",mumps->id.ICNTL(1));CHKERRQ(ierr);
1178a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg): %d \n",mumps->id.ICNTL(2));CHKERRQ(ierr);
1179a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):   %d \n",mumps->id.ICNTL(3));CHKERRQ(ierr);
1180a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):        %d \n",mumps->id.ICNTL(4));CHKERRQ(ierr);
1181a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):         %d \n",mumps->id.ICNTL(5));CHKERRQ(ierr);
1182a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):        %d \n",mumps->id.ICNTL(6));CHKERRQ(ierr);
1183a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (sequentia matrix ordering):%d \n",mumps->id.ICNTL(7));CHKERRQ(ierr);
1184a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scalling strategy):        %d \n",mumps->id.ICNTL(8));CHKERRQ(ierr);
1185a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements):  %d \n",mumps->id.ICNTL(10));CHKERRQ(ierr);
1186a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):          %d \n",mumps->id.ICNTL(11));CHKERRQ(ierr);
1187a5e57a09SHong Zhang       if (mumps->id.ICNTL(11)>0) {
1188a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(4) (inf norm of input mat):        %g\n",mumps->id.RINFOG(4));CHKERRQ(ierr);
1189a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(5) (inf norm of solution):         %g\n",mumps->id.RINFOG(5));CHKERRQ(ierr);
1190a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(6) (inf norm of residual):         %g\n",mumps->id.RINFOG(6));CHKERRQ(ierr);
1191a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8));CHKERRQ(ierr);
1192a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(9) (error estimate):               %g \n",mumps->id.RINFOG(9));CHKERRQ(ierr);
1193a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11));CHKERRQ(ierr);
1194f6c57405SHong Zhang       }
1195a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",mumps->id.ICNTL(12));CHKERRQ(ierr);
1196a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",mumps->id.ICNTL(13));CHKERRQ(ierr);
1197a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",mumps->id.ICNTL(14));CHKERRQ(ierr);
1198f6c57405SHong Zhang       /* ICNTL(15-17) not used */
1199a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",mumps->id.ICNTL(18));CHKERRQ(ierr);
1200a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Shur complement info):                       %d \n",mumps->id.ICNTL(19));CHKERRQ(ierr);
1201a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",mumps->id.ICNTL(20));CHKERRQ(ierr);
1202a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (somumpstion struct):                            %d \n",mumps->id.ICNTL(21));CHKERRQ(ierr);
1203a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",mumps->id.ICNTL(22));CHKERRQ(ierr);
1204a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",mumps->id.ICNTL(23));CHKERRQ(ierr);
1205c0165424SHong Zhang 
1206a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",mumps->id.ICNTL(24));CHKERRQ(ierr);
1207a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",mumps->id.ICNTL(25));CHKERRQ(ierr);
1208a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",mumps->id.ICNTL(26));CHKERRQ(ierr);
1209a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",mumps->id.ICNTL(27));CHKERRQ(ierr);
1210a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(28) (use parallel or sequential ordering):        %d \n",mumps->id.ICNTL(28));CHKERRQ(ierr);
1211a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(29) (parallel ordering):                          %d \n",mumps->id.ICNTL(29));CHKERRQ(ierr);
121242179a6aSHong Zhang 
1213a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(30) (user-specified set of entries in inv(A)):    %d \n",mumps->id.ICNTL(30));CHKERRQ(ierr);
1214a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(31) (factors is discarded in the solve phase):    %d \n",mumps->id.ICNTL(31));CHKERRQ(ierr);
1215a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(33) (compute determinant):                        %d \n",mumps->id.ICNTL(33));CHKERRQ(ierr);
1216f6c57405SHong Zhang 
1217a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",mumps->id.CNTL(1));CHKERRQ(ierr);
1218a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2));CHKERRQ(ierr);
1219a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absomumpste pivoting threshold):      %g \n",mumps->id.CNTL(3));CHKERRQ(ierr);
1220a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(4) (vamumpse of static pivoting):         %g \n",mumps->id.CNTL(4));CHKERRQ(ierr);
1221a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",mumps->id.CNTL(5));CHKERRQ(ierr);
1222f6c57405SHong Zhang 
1223f6c57405SHong Zhang       /* infomation local to each processor */
122434ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);
12257b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
1226a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %g \n",mumps->myid,mumps->id.RINFO(1));CHKERRQ(ierr);
122734ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
122834ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);
1229a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(2));CHKERRQ(ierr);
123034ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
123134ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);
1232a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(3));CHKERRQ(ierr);
123334ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1234f6c57405SHong Zhang 
123534ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr);
1236a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] %d \n",mumps->myid,mumps->id.INFO(15));CHKERRQ(ierr);
123734ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1238f6c57405SHong Zhang 
123934ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr);
1240a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(16));CHKERRQ(ierr);
124134ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1242f6c57405SHong Zhang 
124334ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);
1244a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(23));CHKERRQ(ierr);
124534ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
12467b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
1247f6c57405SHong Zhang 
1248a5e57a09SHong Zhang       if (!mumps->myid) { /* information from the host */
1249a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",mumps->id.RINFOG(1));CHKERRQ(ierr);
1250a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",mumps->id.RINFOG(2));CHKERRQ(ierr);
1251a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",mumps->id.RINFOG(3));CHKERRQ(ierr);
1252a5e57a09SHong 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);
1253f6c57405SHong Zhang 
1254a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(3));CHKERRQ(ierr);
1255a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(4));CHKERRQ(ierr);
1256a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",mumps->id.INFOG(5));CHKERRQ(ierr);
1257a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",mumps->id.INFOG(6));CHKERRQ(ierr);
1258a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively use after analysis): %d \n",mumps->id.INFOG(7));CHKERRQ(ierr);
1259a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",mumps->id.INFOG(8));CHKERRQ(ierr);
1260a5e57a09SHong 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);
1261a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",mumps->id.INFOG(10));CHKERRQ(ierr);
1262a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",mumps->id.INFOG(11));CHKERRQ(ierr);
1263a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",mumps->id.INFOG(12));CHKERRQ(ierr);
1264a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",mumps->id.INFOG(13));CHKERRQ(ierr);
1265a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",mumps->id.INFOG(14));CHKERRQ(ierr);
1266a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",mumps->id.INFOG(15));CHKERRQ(ierr);
1267a5e57a09SHong 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);
1268a5e57a09SHong 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);
1269a5e57a09SHong 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);
1270a5e57a09SHong 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);
1271a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",mumps->id.INFOG(20));CHKERRQ(ierr);
1272a5e57a09SHong 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);
1273a5e57a09SHong 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);
1274a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",mumps->id.INFOG(23));CHKERRQ(ierr);
1275a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",mumps->id.INFOG(24));CHKERRQ(ierr);
1276a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",mumps->id.INFOG(25));CHKERRQ(ierr);
127740d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(28) (after factorization: number of null pivots encountered): %d\n",mumps->id.INFOG(28));CHKERRQ(ierr);
127840d435e3SHong 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);
127940d435e3SHong 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);
128040d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(32) (after analysis: type of analysis done): %d\n",mumps->id.INFOG(32));CHKERRQ(ierr);
128140d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(33) (value used for ICNTL(8)): %d\n",mumps->id.INFOG(33));CHKERRQ(ierr);
128240d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(34) (exponent of the determinant if determinant is requested): %d\n",mumps->id.INFOG(34));CHKERRQ(ierr);
1283f6c57405SHong Zhang       }
1284f6c57405SHong Zhang     }
1285cb828f0fSHong Zhang   }
1286f6c57405SHong Zhang   PetscFunctionReturn(0);
1287f6c57405SHong Zhang }
1288f6c57405SHong Zhang 
128935bd34faSBarry Smith #undef __FUNCT__
129035bd34faSBarry Smith #define __FUNCT__ "MatGetInfo_MUMPS"
129135bd34faSBarry Smith PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
129235bd34faSBarry Smith {
1293cb828f0fSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)A->spptr;
129435bd34faSBarry Smith 
129535bd34faSBarry Smith   PetscFunctionBegin;
129635bd34faSBarry Smith   info->block_size        = 1.0;
1297cb828f0fSHong Zhang   info->nz_allocated      = mumps->id.INFOG(20);
1298cb828f0fSHong Zhang   info->nz_used           = mumps->id.INFOG(20);
129935bd34faSBarry Smith   info->nz_unneeded       = 0.0;
130035bd34faSBarry Smith   info->assemblies        = 0.0;
130135bd34faSBarry Smith   info->mallocs           = 0.0;
130235bd34faSBarry Smith   info->memory            = 0.0;
130335bd34faSBarry Smith   info->fill_ratio_given  = 0;
130435bd34faSBarry Smith   info->fill_ratio_needed = 0;
130535bd34faSBarry Smith   info->factor_mallocs    = 0;
130635bd34faSBarry Smith   PetscFunctionReturn(0);
130735bd34faSBarry Smith }
130835bd34faSBarry Smith 
13095ccb76cbSHong Zhang /* -------------------------------------------------------------------------------------------*/
13105ccb76cbSHong Zhang #undef __FUNCT__
13115ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl_MUMPS"
13125ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival)
13135ccb76cbSHong Zhang {
1314a5e57a09SHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
13155ccb76cbSHong Zhang 
13165ccb76cbSHong Zhang   PetscFunctionBegin;
1317a5e57a09SHong Zhang   mumps->id.ICNTL(icntl) = ival;
13185ccb76cbSHong Zhang   PetscFunctionReturn(0);
13195ccb76cbSHong Zhang }
13205ccb76cbSHong Zhang 
13215ccb76cbSHong Zhang #undef __FUNCT__
1322bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetIcntl_MUMPS"
1323bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt *ival)
1324bc6112feSHong Zhang {
1325bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1326bc6112feSHong Zhang 
1327bc6112feSHong Zhang   PetscFunctionBegin;
1328bc6112feSHong Zhang   *ival = mumps->id.ICNTL(icntl);
1329bc6112feSHong Zhang   PetscFunctionReturn(0);
1330bc6112feSHong Zhang }
1331bc6112feSHong Zhang 
1332bc6112feSHong Zhang #undef __FUNCT__
13335ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl"
13345ccb76cbSHong Zhang /*@
13355ccb76cbSHong Zhang   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()
13365ccb76cbSHong Zhang 
13375ccb76cbSHong Zhang    Logically Collective on Mat
13385ccb76cbSHong Zhang 
13395ccb76cbSHong Zhang    Input Parameters:
13405ccb76cbSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
13415ccb76cbSHong Zhang .  icntl - index of MUMPS parameter array ICNTL()
13425ccb76cbSHong Zhang -  ival - value of MUMPS ICNTL(icntl)
13435ccb76cbSHong Zhang 
13445ccb76cbSHong Zhang   Options Database:
13455ccb76cbSHong Zhang .   -mat_mumps_icntl_<icntl> <ival>
13465ccb76cbSHong Zhang 
13475ccb76cbSHong Zhang    Level: beginner
13485ccb76cbSHong Zhang 
13495ccb76cbSHong Zhang    References: MUMPS Users' Guide
13505ccb76cbSHong Zhang 
13515ccb76cbSHong Zhang .seealso: MatGetFactor()
13525ccb76cbSHong Zhang @*/
13535ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival)
13545ccb76cbSHong Zhang {
13555ccb76cbSHong Zhang   PetscErrorCode ierr;
13565ccb76cbSHong Zhang 
13575ccb76cbSHong Zhang   PetscFunctionBegin;
13585ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
13595ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,ival,3);
13605ccb76cbSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr);
13615ccb76cbSHong Zhang   PetscFunctionReturn(0);
13625ccb76cbSHong Zhang }
13635ccb76cbSHong Zhang 
1364bc6112feSHong Zhang #undef __FUNCT__
1365bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetIcntl"
1366a21f80fcSHong Zhang /*@
1367a21f80fcSHong Zhang   MatMumpsGetIcntl - Get MUMPS parameter ICNTL()
1368a21f80fcSHong Zhang 
1369a21f80fcSHong Zhang    Logically Collective on Mat
1370a21f80fcSHong Zhang 
1371a21f80fcSHong Zhang    Input Parameters:
1372a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1373a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array ICNTL()
1374a21f80fcSHong Zhang 
1375a21f80fcSHong Zhang   Output Parameter:
1376a21f80fcSHong Zhang .  ival - value of MUMPS ICNTL(icntl)
1377a21f80fcSHong Zhang 
1378a21f80fcSHong Zhang    Level: beginner
1379a21f80fcSHong Zhang 
1380a21f80fcSHong Zhang    References: MUMPS Users' Guide
1381a21f80fcSHong Zhang 
1382a21f80fcSHong Zhang .seealso: MatGetFactor()
1383a21f80fcSHong Zhang @*/
1384bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl(Mat F,PetscInt icntl,PetscInt *ival)
1385bc6112feSHong Zhang {
1386bc6112feSHong Zhang   PetscErrorCode ierr;
1387bc6112feSHong Zhang 
1388bc6112feSHong Zhang   PetscFunctionBegin;
1389bc6112feSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
1390bc6112feSHong Zhang   PetscValidIntPointer(ival,3);
1391bc6112feSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetIcntl_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
1392bc6112feSHong Zhang   PetscFunctionReturn(0);
1393bc6112feSHong Zhang }
1394bc6112feSHong Zhang 
13958928b65cSHong Zhang /* -------------------------------------------------------------------------------------------*/
13968928b65cSHong Zhang #undef __FUNCT__
13978928b65cSHong Zhang #define __FUNCT__ "MatMumpsSetCntl_MUMPS"
13988928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal val)
13998928b65cSHong Zhang {
14008928b65cSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
14018928b65cSHong Zhang 
14028928b65cSHong Zhang   PetscFunctionBegin;
14038928b65cSHong Zhang   mumps->id.CNTL(icntl) = val;
14048928b65cSHong Zhang   PetscFunctionReturn(0);
14058928b65cSHong Zhang }
14068928b65cSHong Zhang 
14078928b65cSHong Zhang #undef __FUNCT__
1408bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetCntl_MUMPS"
1409bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal *val)
1410bc6112feSHong Zhang {
1411bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1412bc6112feSHong Zhang 
1413bc6112feSHong Zhang   PetscFunctionBegin;
1414bc6112feSHong Zhang   *val = mumps->id.CNTL(icntl);
1415bc6112feSHong Zhang   PetscFunctionReturn(0);
1416bc6112feSHong Zhang }
1417bc6112feSHong Zhang 
1418bc6112feSHong Zhang #undef __FUNCT__
14198928b65cSHong Zhang #define __FUNCT__ "MatMumpsSetCntl"
14208928b65cSHong Zhang /*@
14218928b65cSHong Zhang   MatMumpsSetCntl - Set MUMPS parameter CNTL()
14228928b65cSHong Zhang 
14238928b65cSHong Zhang    Logically Collective on Mat
14248928b65cSHong Zhang 
14258928b65cSHong Zhang    Input Parameters:
14268928b65cSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
14278928b65cSHong Zhang .  icntl - index of MUMPS parameter array CNTL()
14288928b65cSHong Zhang -  val - value of MUMPS CNTL(icntl)
14298928b65cSHong Zhang 
14308928b65cSHong Zhang   Options Database:
14318928b65cSHong Zhang .   -mat_mumps_cntl_<icntl> <val>
14328928b65cSHong Zhang 
14338928b65cSHong Zhang    Level: beginner
14348928b65cSHong Zhang 
14358928b65cSHong Zhang    References: MUMPS Users' Guide
14368928b65cSHong Zhang 
14378928b65cSHong Zhang .seealso: MatGetFactor()
14388928b65cSHong Zhang @*/
14398928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl(Mat F,PetscInt icntl,PetscReal val)
14408928b65cSHong Zhang {
14418928b65cSHong Zhang   PetscErrorCode ierr;
14428928b65cSHong Zhang 
14438928b65cSHong Zhang   PetscFunctionBegin;
14448928b65cSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
1445bc6112feSHong Zhang   PetscValidLogicalCollectiveReal(F,val,3);
14468928b65cSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsSetCntl_C",(Mat,PetscInt,PetscReal),(F,icntl,val));CHKERRQ(ierr);
14478928b65cSHong Zhang   PetscFunctionReturn(0);
14488928b65cSHong Zhang }
14498928b65cSHong Zhang 
1450bc6112feSHong Zhang #undef __FUNCT__
1451bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetCntl"
1452a21f80fcSHong Zhang /*@
1453a21f80fcSHong Zhang   MatMumpsGetCntl - Get MUMPS parameter CNTL()
1454a21f80fcSHong Zhang 
1455a21f80fcSHong Zhang    Logically Collective on Mat
1456a21f80fcSHong Zhang 
1457a21f80fcSHong Zhang    Input Parameters:
1458a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1459a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array CNTL()
1460a21f80fcSHong Zhang 
1461a21f80fcSHong Zhang   Output Parameter:
1462a21f80fcSHong Zhang .  val - value of MUMPS CNTL(icntl)
1463a21f80fcSHong Zhang 
1464a21f80fcSHong Zhang    Level: beginner
1465a21f80fcSHong Zhang 
1466a21f80fcSHong Zhang    References: MUMPS Users' Guide
1467a21f80fcSHong Zhang 
1468a21f80fcSHong Zhang .seealso: MatGetFactor()
1469a21f80fcSHong Zhang @*/
1470bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl(Mat F,PetscInt icntl,PetscReal *val)
1471bc6112feSHong Zhang {
1472bc6112feSHong Zhang   PetscErrorCode ierr;
1473bc6112feSHong Zhang 
1474bc6112feSHong Zhang   PetscFunctionBegin;
1475bc6112feSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
1476bc6112feSHong Zhang   PetscValidRealPointer(val,3);
1477bc6112feSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetCntl_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
1478bc6112feSHong Zhang   PetscFunctionReturn(0);
1479bc6112feSHong Zhang }
1480bc6112feSHong Zhang 
1481bc6112feSHong Zhang #undef __FUNCT__
1482ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfo_MUMPS"
1483ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F,PetscInt icntl,PetscInt *info)
1484bc6112feSHong Zhang {
1485bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1486bc6112feSHong Zhang 
1487bc6112feSHong Zhang   PetscFunctionBegin;
1488bc6112feSHong Zhang   *info = mumps->id.INFO(icntl);
1489bc6112feSHong Zhang   PetscFunctionReturn(0);
1490bc6112feSHong Zhang }
1491bc6112feSHong Zhang 
1492bc6112feSHong Zhang #undef __FUNCT__
1493ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfog_MUMPS"
1494ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F,PetscInt icntl,PetscInt *infog)
1495bc6112feSHong Zhang {
1496bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1497bc6112feSHong Zhang 
1498bc6112feSHong Zhang   PetscFunctionBegin;
1499bc6112feSHong Zhang   *infog = mumps->id.INFOG(icntl);
1500bc6112feSHong Zhang   PetscFunctionReturn(0);
1501bc6112feSHong Zhang }
1502bc6112feSHong Zhang 
1503bc6112feSHong Zhang #undef __FUNCT__
1504ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfo_MUMPS"
1505ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfo)
1506bc6112feSHong Zhang {
1507bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1508bc6112feSHong Zhang 
1509bc6112feSHong Zhang   PetscFunctionBegin;
1510bc6112feSHong Zhang   *rinfo = mumps->id.RINFO(icntl);
1511bc6112feSHong Zhang   PetscFunctionReturn(0);
1512bc6112feSHong Zhang }
1513bc6112feSHong Zhang 
1514bc6112feSHong Zhang #undef __FUNCT__
1515ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfog_MUMPS"
1516ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfog)
1517bc6112feSHong Zhang {
1518bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1519bc6112feSHong Zhang 
1520bc6112feSHong Zhang   PetscFunctionBegin;
1521bc6112feSHong Zhang   *rinfog = mumps->id.RINFOG(icntl);
1522bc6112feSHong Zhang   PetscFunctionReturn(0);
1523bc6112feSHong Zhang }
1524bc6112feSHong Zhang 
1525bc6112feSHong Zhang #undef __FUNCT__
1526ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfo"
1527a21f80fcSHong Zhang /*@
1528a21f80fcSHong Zhang   MatMumpsGetInfo - Get MUMPS parameter INFO()
1529a21f80fcSHong Zhang 
1530a21f80fcSHong Zhang    Logically Collective on Mat
1531a21f80fcSHong Zhang 
1532a21f80fcSHong Zhang    Input Parameters:
1533a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1534a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array INFO()
1535a21f80fcSHong Zhang 
1536a21f80fcSHong Zhang   Output Parameter:
1537a21f80fcSHong Zhang .  ival - value of MUMPS INFO(icntl)
1538a21f80fcSHong Zhang 
1539a21f80fcSHong Zhang    Level: beginner
1540a21f80fcSHong Zhang 
1541a21f80fcSHong Zhang    References: MUMPS Users' Guide
1542a21f80fcSHong Zhang 
1543a21f80fcSHong Zhang .seealso: MatGetFactor()
1544a21f80fcSHong Zhang @*/
1545ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo(Mat F,PetscInt icntl,PetscInt *ival)
1546bc6112feSHong Zhang {
1547bc6112feSHong Zhang   PetscErrorCode ierr;
1548bc6112feSHong Zhang 
1549bc6112feSHong Zhang   PetscFunctionBegin;
1550ca810319SHong Zhang   PetscValidIntPointer(ival,3);
1551ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetInfo_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
1552bc6112feSHong Zhang   PetscFunctionReturn(0);
1553bc6112feSHong Zhang }
1554bc6112feSHong Zhang 
1555bc6112feSHong Zhang #undef __FUNCT__
1556ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfog"
1557a21f80fcSHong Zhang /*@
1558a21f80fcSHong Zhang   MatMumpsGetInfog - Get MUMPS parameter INFOG()
1559a21f80fcSHong Zhang 
1560a21f80fcSHong Zhang    Logically Collective on Mat
1561a21f80fcSHong Zhang 
1562a21f80fcSHong Zhang    Input Parameters:
1563a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1564a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array INFOG()
1565a21f80fcSHong Zhang 
1566a21f80fcSHong Zhang   Output Parameter:
1567a21f80fcSHong Zhang .  ival - value of MUMPS INFOG(icntl)
1568a21f80fcSHong Zhang 
1569a21f80fcSHong Zhang    Level: beginner
1570a21f80fcSHong Zhang 
1571a21f80fcSHong Zhang    References: MUMPS Users' Guide
1572a21f80fcSHong Zhang 
1573a21f80fcSHong Zhang .seealso: MatGetFactor()
1574a21f80fcSHong Zhang @*/
1575ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog(Mat F,PetscInt icntl,PetscInt *ival)
1576bc6112feSHong Zhang {
1577bc6112feSHong Zhang   PetscErrorCode ierr;
1578bc6112feSHong Zhang 
1579bc6112feSHong Zhang   PetscFunctionBegin;
1580ca810319SHong Zhang   PetscValidIntPointer(ival,3);
1581ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetInfog_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
1582bc6112feSHong Zhang   PetscFunctionReturn(0);
1583bc6112feSHong Zhang }
1584bc6112feSHong Zhang 
1585bc6112feSHong Zhang #undef __FUNCT__
1586ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfo"
1587a21f80fcSHong Zhang /*@
1588a21f80fcSHong Zhang   MatMumpsGetRinfo - Get MUMPS parameter RINFO()
1589a21f80fcSHong Zhang 
1590a21f80fcSHong Zhang    Logically Collective on Mat
1591a21f80fcSHong Zhang 
1592a21f80fcSHong Zhang    Input Parameters:
1593a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1594a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array RINFO()
1595a21f80fcSHong Zhang 
1596a21f80fcSHong Zhang   Output Parameter:
1597a21f80fcSHong Zhang .  val - value of MUMPS RINFO(icntl)
1598a21f80fcSHong Zhang 
1599a21f80fcSHong Zhang    Level: beginner
1600a21f80fcSHong Zhang 
1601a21f80fcSHong Zhang    References: MUMPS Users' Guide
1602a21f80fcSHong Zhang 
1603a21f80fcSHong Zhang .seealso: MatGetFactor()
1604a21f80fcSHong Zhang @*/
1605ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo(Mat F,PetscInt icntl,PetscReal *val)
1606bc6112feSHong Zhang {
1607bc6112feSHong Zhang   PetscErrorCode ierr;
1608bc6112feSHong Zhang 
1609bc6112feSHong Zhang   PetscFunctionBegin;
1610bc6112feSHong Zhang   PetscValidRealPointer(val,3);
1611ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetRinfo_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
1612bc6112feSHong Zhang   PetscFunctionReturn(0);
1613bc6112feSHong Zhang }
1614bc6112feSHong Zhang 
1615bc6112feSHong Zhang #undef __FUNCT__
1616ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfog"
1617a21f80fcSHong Zhang /*@
1618a21f80fcSHong Zhang   MatMumpsGetRinfog - Get MUMPS parameter RINFOG()
1619a21f80fcSHong Zhang 
1620a21f80fcSHong Zhang    Logically Collective on Mat
1621a21f80fcSHong Zhang 
1622a21f80fcSHong Zhang    Input Parameters:
1623a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1624a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array RINFOG()
1625a21f80fcSHong Zhang 
1626a21f80fcSHong Zhang   Output Parameter:
1627a21f80fcSHong Zhang .  val - value of MUMPS RINFOG(icntl)
1628a21f80fcSHong Zhang 
1629a21f80fcSHong Zhang    Level: beginner
1630a21f80fcSHong Zhang 
1631a21f80fcSHong Zhang    References: MUMPS Users' Guide
1632a21f80fcSHong Zhang 
1633a21f80fcSHong Zhang .seealso: MatGetFactor()
1634a21f80fcSHong Zhang @*/
1635ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog(Mat F,PetscInt icntl,PetscReal *val)
1636bc6112feSHong Zhang {
1637bc6112feSHong Zhang   PetscErrorCode ierr;
1638bc6112feSHong Zhang 
1639bc6112feSHong Zhang   PetscFunctionBegin;
1640bc6112feSHong Zhang   PetscValidRealPointer(val,3);
1641ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetRinfog_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
1642bc6112feSHong Zhang   PetscFunctionReturn(0);
1643bc6112feSHong Zhang }
1644bc6112feSHong Zhang 
164524b6179bSKris Buschelman /*MC
16462692d6eeSBarry Smith   MATSOLVERMUMPS -  A matrix type providing direct solvers (LU and Cholesky) for
164724b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
164824b6179bSKris Buschelman 
164941c8de11SBarry Smith   Works with MATAIJ and MATSBAIJ matrices
165024b6179bSKris Buschelman 
165124b6179bSKris Buschelman   Options Database Keys:
1652fb8376fbSHong Zhang + -mat_mumps_icntl_4 <0,...,4> - print level
165324b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
165464e6c443SBarry Smith . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guidec)
165524b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
165624b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
165794b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
165824b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
165924b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
166024b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
166124b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
166224b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
166324b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
166424b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
166524b6179bSKris Buschelman 
166624b6179bSKris Buschelman   Level: beginner
166724b6179bSKris Buschelman 
166841c8de11SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage
166941c8de11SBarry Smith 
167024b6179bSKris Buschelman M*/
167124b6179bSKris Buschelman 
167235bd34faSBarry Smith #undef __FUNCT__
167335bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_mumps"
1674f7a08781SBarry Smith static PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type)
167535bd34faSBarry Smith {
167635bd34faSBarry Smith   PetscFunctionBegin;
16772692d6eeSBarry Smith   *type = MATSOLVERMUMPS;
167835bd34faSBarry Smith   PetscFunctionReturn(0);
167935bd34faSBarry Smith }
168035bd34faSBarry Smith 
1681bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI AIJ matrices */
16822877fffaSHong Zhang #undef __FUNCT__
1683bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_aij_mumps"
16848cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
16852877fffaSHong Zhang {
16862877fffaSHong Zhang   Mat            B;
16872877fffaSHong Zhang   PetscErrorCode ierr;
16882877fffaSHong Zhang   Mat_MUMPS      *mumps;
1689ace3abfcSBarry Smith   PetscBool      isSeqAIJ;
16902877fffaSHong Zhang 
16912877fffaSHong Zhang   PetscFunctionBegin;
16922877fffaSHong Zhang   /* Create the factorization matrix */
1693251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr);
1694ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
16952877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
16962877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1697bccb9932SShri Abhyankar   if (isSeqAIJ) {
16980298fd71SBarry Smith     ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr);
1699bccb9932SShri Abhyankar   } else {
17000298fd71SBarry Smith     ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr);
1701bccb9932SShri Abhyankar   }
17022877fffaSHong Zhang 
1703b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
17042205254eSKarl Rupp 
17052877fffaSHong Zhang   B->ops->view        = MatView_MUMPS;
170635bd34faSBarry Smith   B->ops->getinfo     = MatGetInfo_MUMPS;
170720be8e61SHong Zhang   B->ops->getdiagonal = MatGetDiagonal_MUMPS;
17082205254eSKarl Rupp 
1709bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
1710bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
1711bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
1712bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
1713bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
1714bc6112feSHong Zhang 
1715ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
1716ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
1717ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
1718ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
1719450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
1720450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
1721d5f3da31SBarry Smith     B->factortype            = MAT_FACTOR_LU;
1722bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
1723bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
1724746480a1SHong Zhang     mumps->sym = 0;
1725dcd589f8SShri Abhyankar   } else {
172667877ebaSShri Abhyankar     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
1727450b117fSShri Abhyankar     B->factortype                  = MAT_FACTOR_CHOLESKY;
1728bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
1729bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
17306fdc2a6dSBarry Smith     if (A->spd_set && A->spd) mumps->sym = 1;
17316fdc2a6dSBarry Smith     else                      mumps->sym = 2;
1732450b117fSShri Abhyankar   }
17332877fffaSHong Zhang 
17342877fffaSHong Zhang   mumps->isAIJ    = PETSC_TRUE;
1735bf0cc555SLisandro Dalcin   mumps->Destroy  = B->ops->destroy;
17362877fffaSHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
17372877fffaSHong Zhang   B->spptr        = (void*)mumps;
17382205254eSKarl Rupp 
1739f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1740746480a1SHong Zhang 
17412877fffaSHong Zhang   *F = B;
17422877fffaSHong Zhang   PetscFunctionReturn(0);
17432877fffaSHong Zhang }
17442877fffaSHong Zhang 
1745bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI SBAIJ matrices */
17462877fffaSHong Zhang #undef __FUNCT__
1747bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_sbaij_mumps"
17488cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
17492877fffaSHong Zhang {
17502877fffaSHong Zhang   Mat            B;
17512877fffaSHong Zhang   PetscErrorCode ierr;
17522877fffaSHong Zhang   Mat_MUMPS      *mumps;
1753ace3abfcSBarry Smith   PetscBool      isSeqSBAIJ;
17542877fffaSHong Zhang 
17552877fffaSHong Zhang   PetscFunctionBegin;
1756ce94432eSBarry Smith   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
1757ce94432eSBarry 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");
1758251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
17592877fffaSHong Zhang   /* Create the factorization matrix */
1760ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
17612877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
17622877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1763b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
1764bccb9932SShri Abhyankar   if (isSeqSBAIJ) {
17650298fd71SBarry Smith     ierr = MatSeqSBAIJSetPreallocation(B,1,0,NULL);CHKERRQ(ierr);
17662205254eSKarl Rupp 
176716ebf90aSShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
1768dcd589f8SShri Abhyankar   } else {
17690298fd71SBarry Smith     ierr = MatMPISBAIJSetPreallocation(B,1,0,NULL,0,NULL);CHKERRQ(ierr);
17702205254eSKarl Rupp 
1771bccb9932SShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
1772bccb9932SShri Abhyankar   }
1773bccb9932SShri Abhyankar 
177467877ebaSShri Abhyankar   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
1775bccb9932SShri Abhyankar   B->ops->view                   = MatView_MUMPS;
177620be8e61SHong Zhang   B->ops->getdiagonal            = MatGetDiagonal_MUMPS;
17772205254eSKarl Rupp 
1778bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
1779b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
1780b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
1781b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
1782b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
1783bc6112feSHong Zhang 
1784ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
1785ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
1786ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
1787ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
17882205254eSKarl Rupp 
1789f4762488SHong Zhang   B->factortype = MAT_FACTOR_CHOLESKY;
17906fdc2a6dSBarry Smith   if (A->spd_set && A->spd) mumps->sym = 1;
17916fdc2a6dSBarry Smith   else                      mumps->sym = 2;
1792a214ac2aSShri Abhyankar 
1793bccb9932SShri Abhyankar   mumps->isAIJ    = PETSC_FALSE;
1794bf0cc555SLisandro Dalcin   mumps->Destroy  = B->ops->destroy;
1795f3c0ef26SHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
17962877fffaSHong Zhang   B->spptr        = (void*)mumps;
17972205254eSKarl Rupp 
1798f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1799746480a1SHong Zhang 
18002877fffaSHong Zhang   *F = B;
18012877fffaSHong Zhang   PetscFunctionReturn(0);
18022877fffaSHong Zhang }
180397969023SHong Zhang 
1804450b117fSShri Abhyankar #undef __FUNCT__
1805bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_baij_mumps"
18068cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
180767877ebaSShri Abhyankar {
180867877ebaSShri Abhyankar   Mat            B;
180967877ebaSShri Abhyankar   PetscErrorCode ierr;
181067877ebaSShri Abhyankar   Mat_MUMPS      *mumps;
1811ace3abfcSBarry Smith   PetscBool      isSeqBAIJ;
181267877ebaSShri Abhyankar 
181367877ebaSShri Abhyankar   PetscFunctionBegin;
181467877ebaSShri Abhyankar   /* Create the factorization matrix */
1815251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr);
1816ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
181767877ebaSShri Abhyankar   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
181867877ebaSShri Abhyankar   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1819bccb9932SShri Abhyankar   if (isSeqBAIJ) {
18200298fd71SBarry Smith     ierr = MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,NULL);CHKERRQ(ierr);
1821bccb9932SShri Abhyankar   } else {
18220298fd71SBarry Smith     ierr = MatMPIBAIJSetPreallocation(B,A->rmap->bs,0,NULL,0,NULL);CHKERRQ(ierr);
1823bccb9932SShri Abhyankar   }
1824450b117fSShri Abhyankar 
1825b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
1826450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
1827450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
1828450b117fSShri Abhyankar     B->factortype            = MAT_FACTOR_LU;
1829bccb9932SShri Abhyankar     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
1830bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
1831746480a1SHong Zhang     mumps->sym = 0;
1832f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n");
1833bccb9932SShri Abhyankar 
1834450b117fSShri Abhyankar   B->ops->view        = MatView_MUMPS;
183520be8e61SHong Zhang   B->ops->getdiagonal = MatGetDiagonal_MUMPS;
18362205254eSKarl Rupp 
1837bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
1838bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
1839bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
1840bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
1841bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
1842bc6112feSHong Zhang 
1843ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
1844ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
1845ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
1846ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
1847450b117fSShri Abhyankar 
1848450b117fSShri Abhyankar   mumps->isAIJ    = PETSC_TRUE;
1849bf0cc555SLisandro Dalcin   mumps->Destroy  = B->ops->destroy;
1850450b117fSShri Abhyankar   B->ops->destroy = MatDestroy_MUMPS;
1851450b117fSShri Abhyankar   B->spptr        = (void*)mumps;
18522205254eSKarl Rupp 
1853f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1854746480a1SHong Zhang 
1855450b117fSShri Abhyankar   *F = B;
1856450b117fSShri Abhyankar   PetscFunctionReturn(0);
1857450b117fSShri Abhyankar }
185842c9c57cSBarry Smith 
185942c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mumps(Mat,MatFactorType,Mat*);
186042c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_baij_mumps(Mat,MatFactorType,Mat*);
186142c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat,MatFactorType,Mat*);
186242c9c57cSBarry Smith 
186342c9c57cSBarry Smith #undef __FUNCT__
186442c9c57cSBarry Smith #define __FUNCT__ "MatSolverPackageRegister_MUMPS"
186529b38603SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_MUMPS(void)
186642c9c57cSBarry Smith {
186742c9c57cSBarry Smith   PetscErrorCode ierr;
186842c9c57cSBarry Smith 
186942c9c57cSBarry Smith   PetscFunctionBegin;
187042c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ,         MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr);
187142c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ,         MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr);
187242c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ,         MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr);
187342c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ,         MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr);
187442c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPISBAIJ,         MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr);
187542c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ,         MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr);
187642c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ,         MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr);
187742c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ,         MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr);
187842c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ,         MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr);
187942c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQSBAIJ,         MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr);
188042c9c57cSBarry Smith   PetscFunctionReturn(0);
188142c9c57cSBarry Smith }
188242c9c57cSBarry Smith 
1883