xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision a21f80fc8d342514f0f71e9bcfcf6939aaa9c57d)
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 
8567877ebaSShri Abhyankar 
8667877ebaSShri Abhyankar /* MatConvertToTriples_A_B */
8767877ebaSShri Abhyankar /*convert Petsc matrix to triples: row[nz], col[nz], val[nz] */
88397b6df1SKris Buschelman /*
89397b6df1SKris Buschelman   input:
9067877ebaSShri Abhyankar     A       - matrix in aij,baij or sbaij (bs=1) format
91397b6df1SKris Buschelman     shift   - 0: C style output triple; 1: Fortran style output triple.
92bccb9932SShri Abhyankar     reuse   - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple
93bccb9932SShri Abhyankar               MAT_REUSE_MATRIX:   only the values in v array are updated
94397b6df1SKris Buschelman   output:
95397b6df1SKris Buschelman     nnz     - dim of r, c, and v (number of local nonzero entries of A)
96397b6df1SKris Buschelman     r, c, v - row and col index, matrix values (matrix triples)
97eb9baa12SBarry Smith 
98eb9baa12SBarry Smith   The returned values r, c, and sometimes v are obtained in a single PetscMalloc(). Then in MatDestroy_MUMPS() it is
99eb9baa12SBarry Smith   freed with PetscFree((mumps->irn);  This is not ideal code, the fact that v is ONLY sometimes part of mumps->irn means
100eb9baa12SBarry Smith   that the PetscMalloc() cannot easily be replaced with a PetscMalloc3().
101eb9baa12SBarry Smith 
102397b6df1SKris Buschelman  */
10316ebf90aSShri Abhyankar 
10416ebf90aSShri Abhyankar #undef __FUNCT__
10516ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqaij"
106bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
107b24902e0SBarry Smith {
108185f6596SHong Zhang   const PetscInt *ai,*aj,*ajj,M=A->rmap->n;
10967877ebaSShri Abhyankar   PetscInt       nz,rnz,i,j;
110dfbe8321SBarry Smith   PetscErrorCode ierr;
111c1490034SHong Zhang   PetscInt       *row,*col;
11216ebf90aSShri Abhyankar   Mat_SeqAIJ     *aa=(Mat_SeqAIJ*)A->data;
113397b6df1SKris Buschelman 
114397b6df1SKris Buschelman   PetscFunctionBegin;
11516ebf90aSShri Abhyankar   *v=aa->a;
116bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
1172205254eSKarl Rupp     nz   = aa->nz;
1182205254eSKarl Rupp     ai   = aa->i;
1192205254eSKarl Rupp     aj   = aa->j;
12016ebf90aSShri Abhyankar     *nnz = nz;
121785e854fSJed Brown     ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr);
122185f6596SHong Zhang     col  = row + nz;
123185f6596SHong Zhang 
12416ebf90aSShri Abhyankar     nz = 0;
12516ebf90aSShri Abhyankar     for (i=0; i<M; i++) {
12616ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
12767877ebaSShri Abhyankar       ajj = aj + ai[i];
12867877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
12967877ebaSShri Abhyankar         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
13016ebf90aSShri Abhyankar       }
13116ebf90aSShri Abhyankar     }
13216ebf90aSShri Abhyankar     *r = row; *c = col;
13316ebf90aSShri Abhyankar   }
13416ebf90aSShri Abhyankar   PetscFunctionReturn(0);
13516ebf90aSShri Abhyankar }
136397b6df1SKris Buschelman 
13716ebf90aSShri Abhyankar #undef __FUNCT__
13867877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqbaij_seqaij"
139bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
14067877ebaSShri Abhyankar {
14167877ebaSShri Abhyankar   Mat_SeqBAIJ    *aa=(Mat_SeqBAIJ*)A->data;
14233d57670SJed Brown   const PetscInt *ai,*aj,*ajj,bs2 = aa->bs2;
14333d57670SJed Brown   PetscInt       bs,M,nz,idx=0,rnz,i,j,k,m;
14467877ebaSShri Abhyankar   PetscErrorCode ierr;
14567877ebaSShri Abhyankar   PetscInt       *row,*col;
14667877ebaSShri Abhyankar 
14767877ebaSShri Abhyankar   PetscFunctionBegin;
14833d57670SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
14933d57670SJed Brown   M = A->rmap->N/bs;
150cf3759fdSShri Abhyankar   *v = aa->a;
151bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
152cf3759fdSShri Abhyankar     ai   = aa->i; aj = aa->j;
15367877ebaSShri Abhyankar     nz   = bs2*aa->nz;
15467877ebaSShri Abhyankar     *nnz = nz;
155785e854fSJed Brown     ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr);
156185f6596SHong Zhang     col  = row + nz;
157185f6596SHong Zhang 
15867877ebaSShri Abhyankar     for (i=0; i<M; i++) {
15967877ebaSShri Abhyankar       ajj = aj + ai[i];
16067877ebaSShri Abhyankar       rnz = ai[i+1] - ai[i];
16167877ebaSShri Abhyankar       for (k=0; k<rnz; k++) {
16267877ebaSShri Abhyankar         for (j=0; j<bs; j++) {
16367877ebaSShri Abhyankar           for (m=0; m<bs; m++) {
16467877ebaSShri Abhyankar             row[idx]   = i*bs + m + shift;
165cf3759fdSShri Abhyankar             col[idx++] = bs*(ajj[k]) + j + shift;
16667877ebaSShri Abhyankar           }
16767877ebaSShri Abhyankar         }
16867877ebaSShri Abhyankar       }
16967877ebaSShri Abhyankar     }
170cf3759fdSShri Abhyankar     *r = row; *c = col;
17167877ebaSShri Abhyankar   }
17267877ebaSShri Abhyankar   PetscFunctionReturn(0);
17367877ebaSShri Abhyankar }
17467877ebaSShri Abhyankar 
17567877ebaSShri Abhyankar #undef __FUNCT__
17616ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqsbaij_seqsbaij"
177bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
17816ebf90aSShri Abhyankar {
17967877ebaSShri Abhyankar   const PetscInt *ai, *aj,*ajj,M=A->rmap->n;
18067877ebaSShri Abhyankar   PetscInt       nz,rnz,i,j;
18116ebf90aSShri Abhyankar   PetscErrorCode ierr;
18216ebf90aSShri Abhyankar   PetscInt       *row,*col;
18316ebf90aSShri Abhyankar   Mat_SeqSBAIJ   *aa=(Mat_SeqSBAIJ*)A->data;
18416ebf90aSShri Abhyankar 
18516ebf90aSShri Abhyankar   PetscFunctionBegin;
186882afa5aSHong Zhang   *v = aa->a;
187bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
1882205254eSKarl Rupp     nz   = aa->nz;
1892205254eSKarl Rupp     ai   = aa->i;
1902205254eSKarl Rupp     aj   = aa->j;
1912205254eSKarl Rupp     *v   = aa->a;
19216ebf90aSShri Abhyankar     *nnz = nz;
193785e854fSJed Brown     ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr);
194185f6596SHong Zhang     col  = row + nz;
195185f6596SHong Zhang 
19616ebf90aSShri Abhyankar     nz = 0;
19716ebf90aSShri Abhyankar     for (i=0; i<M; i++) {
19816ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
19967877ebaSShri Abhyankar       ajj = aj + ai[i];
20067877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
20167877ebaSShri Abhyankar         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
20216ebf90aSShri Abhyankar       }
20316ebf90aSShri Abhyankar     }
20416ebf90aSShri Abhyankar     *r = row; *c = col;
20516ebf90aSShri Abhyankar   }
20616ebf90aSShri Abhyankar   PetscFunctionReturn(0);
20716ebf90aSShri Abhyankar }
20816ebf90aSShri Abhyankar 
20916ebf90aSShri Abhyankar #undef __FUNCT__
21016ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqsbaij"
211bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
21216ebf90aSShri Abhyankar {
21367877ebaSShri Abhyankar   const PetscInt    *ai,*aj,*ajj,*adiag,M=A->rmap->n;
21467877ebaSShri Abhyankar   PetscInt          nz,rnz,i,j;
21567877ebaSShri Abhyankar   const PetscScalar *av,*v1;
21616ebf90aSShri Abhyankar   PetscScalar       *val;
21716ebf90aSShri Abhyankar   PetscErrorCode    ierr;
21816ebf90aSShri Abhyankar   PetscInt          *row,*col;
21916ebf90aSShri Abhyankar   Mat_SeqSBAIJ      *aa=(Mat_SeqSBAIJ*)A->data;
22016ebf90aSShri Abhyankar 
22116ebf90aSShri Abhyankar   PetscFunctionBegin;
22216ebf90aSShri Abhyankar   ai   =aa->i; aj=aa->j;av=aa->a;
22316ebf90aSShri Abhyankar   adiag=aa->diag;
224bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
22516ebf90aSShri Abhyankar     nz   = M + (aa->nz-M)/2;
22616ebf90aSShri Abhyankar     *nnz = nz;
227185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
228185f6596SHong Zhang     col  = row + nz;
229185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
230185f6596SHong Zhang 
23116ebf90aSShri Abhyankar     nz = 0;
23216ebf90aSShri Abhyankar     for (i=0; i<M; i++) {
23316ebf90aSShri Abhyankar       rnz = ai[i+1] - adiag[i];
23467877ebaSShri Abhyankar       ajj = aj + adiag[i];
235cf3759fdSShri Abhyankar       v1  = av + adiag[i];
23667877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
23767877ebaSShri Abhyankar         row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
23816ebf90aSShri Abhyankar       }
23916ebf90aSShri Abhyankar     }
24016ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
241397b6df1SKris Buschelman   } else {
24216ebf90aSShri Abhyankar     nz = 0; val = *v;
24316ebf90aSShri Abhyankar     for (i=0; i <M; i++) {
24416ebf90aSShri Abhyankar       rnz = ai[i+1] - adiag[i];
24567877ebaSShri Abhyankar       ajj = aj + adiag[i];
24667877ebaSShri Abhyankar       v1  = av + adiag[i];
24767877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
24867877ebaSShri Abhyankar         val[nz++] = v1[j];
24916ebf90aSShri Abhyankar       }
25016ebf90aSShri Abhyankar     }
25116ebf90aSShri Abhyankar   }
25216ebf90aSShri Abhyankar   PetscFunctionReturn(0);
25316ebf90aSShri Abhyankar }
25416ebf90aSShri Abhyankar 
25516ebf90aSShri Abhyankar #undef __FUNCT__
25616ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpisbaij_mpisbaij"
257bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
25816ebf90aSShri Abhyankar {
25916ebf90aSShri Abhyankar   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
26016ebf90aSShri Abhyankar   PetscErrorCode    ierr;
26116ebf90aSShri Abhyankar   PetscInt          rstart,nz,i,j,jj,irow,countA,countB;
26216ebf90aSShri Abhyankar   PetscInt          *row,*col;
26316ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
26416ebf90aSShri Abhyankar   PetscScalar       *val;
265397b6df1SKris Buschelman   Mat_MPISBAIJ      *mat = (Mat_MPISBAIJ*)A->data;
266397b6df1SKris Buschelman   Mat_SeqSBAIJ      *aa  = (Mat_SeqSBAIJ*)(mat->A)->data;
267397b6df1SKris Buschelman   Mat_SeqBAIJ       *bb  = (Mat_SeqBAIJ*)(mat->B)->data;
26816ebf90aSShri Abhyankar 
26916ebf90aSShri Abhyankar   PetscFunctionBegin;
270d0f46423SBarry Smith   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
271397b6df1SKris Buschelman   av=aa->a; bv=bb->a;
272397b6df1SKris Buschelman 
2732205254eSKarl Rupp   garray = mat->garray;
2742205254eSKarl Rupp 
275bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
27616ebf90aSShri Abhyankar     nz   = aa->nz + bb->nz;
27716ebf90aSShri Abhyankar     *nnz = nz;
278185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
279185f6596SHong Zhang     col  = row + nz;
280185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
281185f6596SHong Zhang 
282397b6df1SKris Buschelman     *r = row; *c = col; *v = val;
283397b6df1SKris Buschelman   } else {
284397b6df1SKris Buschelman     row = *r; col = *c; val = *v;
285397b6df1SKris Buschelman   }
286397b6df1SKris Buschelman 
287028e57e8SHong Zhang   jj = 0; irow = rstart;
288397b6df1SKris Buschelman   for (i=0; i<m; i++) {
289397b6df1SKris Buschelman     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
290397b6df1SKris Buschelman     countA = ai[i+1] - ai[i];
291397b6df1SKris Buschelman     countB = bi[i+1] - bi[i];
292397b6df1SKris Buschelman     bjj    = bj + bi[i];
29316ebf90aSShri Abhyankar     v1     = av + ai[i];
29416ebf90aSShri Abhyankar     v2     = bv + bi[i];
295397b6df1SKris Buschelman 
296397b6df1SKris Buschelman     /* A-part */
297397b6df1SKris Buschelman     for (j=0; j<countA; j++) {
298bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
299397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
300397b6df1SKris Buschelman       }
30116ebf90aSShri Abhyankar       val[jj++] = v1[j];
302397b6df1SKris Buschelman     }
30316ebf90aSShri Abhyankar 
30416ebf90aSShri Abhyankar     /* B-part */
30516ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
306bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
307397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
308397b6df1SKris Buschelman       }
30916ebf90aSShri Abhyankar       val[jj++] = v2[j];
31016ebf90aSShri Abhyankar     }
31116ebf90aSShri Abhyankar     irow++;
31216ebf90aSShri Abhyankar   }
31316ebf90aSShri Abhyankar   PetscFunctionReturn(0);
31416ebf90aSShri Abhyankar }
31516ebf90aSShri Abhyankar 
31616ebf90aSShri Abhyankar #undef __FUNCT__
31716ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpiaij"
318bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
31916ebf90aSShri Abhyankar {
32016ebf90aSShri Abhyankar   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
32116ebf90aSShri Abhyankar   PetscErrorCode    ierr;
32216ebf90aSShri Abhyankar   PetscInt          rstart,nz,i,j,jj,irow,countA,countB;
32316ebf90aSShri Abhyankar   PetscInt          *row,*col;
32416ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
32516ebf90aSShri Abhyankar   PetscScalar       *val;
32616ebf90aSShri Abhyankar   Mat_MPIAIJ        *mat = (Mat_MPIAIJ*)A->data;
32716ebf90aSShri Abhyankar   Mat_SeqAIJ        *aa  = (Mat_SeqAIJ*)(mat->A)->data;
32816ebf90aSShri Abhyankar   Mat_SeqAIJ        *bb  = (Mat_SeqAIJ*)(mat->B)->data;
32916ebf90aSShri Abhyankar 
33016ebf90aSShri Abhyankar   PetscFunctionBegin;
33116ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
33216ebf90aSShri Abhyankar   av=aa->a; bv=bb->a;
33316ebf90aSShri Abhyankar 
3342205254eSKarl Rupp   garray = mat->garray;
3352205254eSKarl Rupp 
336bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
33716ebf90aSShri Abhyankar     nz   = aa->nz + bb->nz;
33816ebf90aSShri Abhyankar     *nnz = nz;
339185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
340185f6596SHong Zhang     col  = row + nz;
341185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
342185f6596SHong Zhang 
34316ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
34416ebf90aSShri Abhyankar   } else {
34516ebf90aSShri Abhyankar     row = *r; col = *c; val = *v;
34616ebf90aSShri Abhyankar   }
34716ebf90aSShri Abhyankar 
34816ebf90aSShri Abhyankar   jj = 0; irow = rstart;
34916ebf90aSShri Abhyankar   for (i=0; i<m; i++) {
35016ebf90aSShri Abhyankar     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
35116ebf90aSShri Abhyankar     countA = ai[i+1] - ai[i];
35216ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
35316ebf90aSShri Abhyankar     bjj    = bj + bi[i];
35416ebf90aSShri Abhyankar     v1     = av + ai[i];
35516ebf90aSShri Abhyankar     v2     = bv + bi[i];
35616ebf90aSShri Abhyankar 
35716ebf90aSShri Abhyankar     /* A-part */
35816ebf90aSShri Abhyankar     for (j=0; j<countA; j++) {
359bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
36016ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
36116ebf90aSShri Abhyankar       }
36216ebf90aSShri Abhyankar       val[jj++] = v1[j];
36316ebf90aSShri Abhyankar     }
36416ebf90aSShri Abhyankar 
36516ebf90aSShri Abhyankar     /* B-part */
36616ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
367bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
36816ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
36916ebf90aSShri Abhyankar       }
37016ebf90aSShri Abhyankar       val[jj++] = v2[j];
37116ebf90aSShri Abhyankar     }
37216ebf90aSShri Abhyankar     irow++;
37316ebf90aSShri Abhyankar   }
37416ebf90aSShri Abhyankar   PetscFunctionReturn(0);
37516ebf90aSShri Abhyankar }
37616ebf90aSShri Abhyankar 
37716ebf90aSShri Abhyankar #undef __FUNCT__
37867877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpibaij_mpiaij"
379bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
38067877ebaSShri Abhyankar {
38167877ebaSShri Abhyankar   Mat_MPIBAIJ       *mat    = (Mat_MPIBAIJ*)A->data;
38267877ebaSShri Abhyankar   Mat_SeqBAIJ       *aa     = (Mat_SeqBAIJ*)(mat->A)->data;
38367877ebaSShri Abhyankar   Mat_SeqBAIJ       *bb     = (Mat_SeqBAIJ*)(mat->B)->data;
38467877ebaSShri Abhyankar   const PetscInt    *ai     = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj;
385d985c460SShri Abhyankar   const PetscInt    *garray = mat->garray,mbs=mat->mbs,rstart=A->rmap->rstart;
38633d57670SJed Brown   const PetscInt    bs2=mat->bs2;
38767877ebaSShri Abhyankar   PetscErrorCode    ierr;
38833d57670SJed Brown   PetscInt          bs,nz,i,j,k,n,jj,irow,countA,countB,idx;
38967877ebaSShri Abhyankar   PetscInt          *row,*col;
39067877ebaSShri Abhyankar   const PetscScalar *av=aa->a, *bv=bb->a,*v1,*v2;
39167877ebaSShri Abhyankar   PetscScalar       *val;
39267877ebaSShri Abhyankar 
39367877ebaSShri Abhyankar   PetscFunctionBegin;
39433d57670SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
395bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
39667877ebaSShri Abhyankar     nz   = bs2*(aa->nz + bb->nz);
39767877ebaSShri Abhyankar     *nnz = nz;
398185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
399185f6596SHong Zhang     col  = row + nz;
400185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
401185f6596SHong Zhang 
40267877ebaSShri Abhyankar     *r = row; *c = col; *v = val;
40367877ebaSShri Abhyankar   } else {
40467877ebaSShri Abhyankar     row = *r; col = *c; val = *v;
40567877ebaSShri Abhyankar   }
40667877ebaSShri Abhyankar 
407d985c460SShri Abhyankar   jj = 0; irow = rstart;
40867877ebaSShri Abhyankar   for (i=0; i<mbs; i++) {
40967877ebaSShri Abhyankar     countA = ai[i+1] - ai[i];
41067877ebaSShri Abhyankar     countB = bi[i+1] - bi[i];
41167877ebaSShri Abhyankar     ajj    = aj + ai[i];
41267877ebaSShri Abhyankar     bjj    = bj + bi[i];
41367877ebaSShri Abhyankar     v1     = av + bs2*ai[i];
41467877ebaSShri Abhyankar     v2     = bv + bs2*bi[i];
41567877ebaSShri Abhyankar 
41667877ebaSShri Abhyankar     idx = 0;
41767877ebaSShri Abhyankar     /* A-part */
41867877ebaSShri Abhyankar     for (k=0; k<countA; k++) {
41967877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
42067877ebaSShri Abhyankar         for (n=0; n<bs; n++) {
421bccb9932SShri Abhyankar           if (reuse == MAT_INITIAL_MATRIX) {
422d985c460SShri Abhyankar             row[jj] = irow + n + shift;
423d985c460SShri Abhyankar             col[jj] = rstart + bs*ajj[k] + j + shift;
42467877ebaSShri Abhyankar           }
42567877ebaSShri Abhyankar           val[jj++] = v1[idx++];
42667877ebaSShri Abhyankar         }
42767877ebaSShri Abhyankar       }
42867877ebaSShri Abhyankar     }
42967877ebaSShri Abhyankar 
43067877ebaSShri Abhyankar     idx = 0;
43167877ebaSShri Abhyankar     /* B-part */
43267877ebaSShri Abhyankar     for (k=0; k<countB; k++) {
43367877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
43467877ebaSShri Abhyankar         for (n=0; n<bs; n++) {
435bccb9932SShri Abhyankar           if (reuse == MAT_INITIAL_MATRIX) {
436d985c460SShri Abhyankar             row[jj] = irow + n + shift;
437d985c460SShri Abhyankar             col[jj] = bs*garray[bjj[k]] + j + shift;
43867877ebaSShri Abhyankar           }
439d985c460SShri Abhyankar           val[jj++] = v2[idx++];
44067877ebaSShri Abhyankar         }
44167877ebaSShri Abhyankar       }
44267877ebaSShri Abhyankar     }
443d985c460SShri Abhyankar     irow += bs;
44467877ebaSShri Abhyankar   }
44567877ebaSShri Abhyankar   PetscFunctionReturn(0);
44667877ebaSShri Abhyankar }
44767877ebaSShri Abhyankar 
44867877ebaSShri Abhyankar #undef __FUNCT__
44916ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpisbaij"
450bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
45116ebf90aSShri Abhyankar {
45216ebf90aSShri Abhyankar   const PetscInt    *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
45316ebf90aSShri Abhyankar   PetscErrorCode    ierr;
454e0bace9bSHong Zhang   PetscInt          rstart,nz,nza,nzb,i,j,jj,irow,countA,countB;
45516ebf90aSShri Abhyankar   PetscInt          *row,*col;
45616ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
45716ebf90aSShri Abhyankar   PetscScalar       *val;
45816ebf90aSShri Abhyankar   Mat_MPIAIJ        *mat =  (Mat_MPIAIJ*)A->data;
45916ebf90aSShri Abhyankar   Mat_SeqAIJ        *aa  =(Mat_SeqAIJ*)(mat->A)->data;
46016ebf90aSShri Abhyankar   Mat_SeqAIJ        *bb  =(Mat_SeqAIJ*)(mat->B)->data;
46116ebf90aSShri Abhyankar 
46216ebf90aSShri Abhyankar   PetscFunctionBegin;
46316ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j; adiag=aa->diag;
46416ebf90aSShri Abhyankar   bi=bb->i; bj=bb->j; garray = mat->garray;
46516ebf90aSShri Abhyankar   av=aa->a; bv=bb->a;
4662205254eSKarl Rupp 
46716ebf90aSShri Abhyankar   rstart = A->rmap->rstart;
46816ebf90aSShri Abhyankar 
469bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
470e0bace9bSHong Zhang     nza = 0;    /* num of upper triangular entries in mat->A, including diagonals */
471e0bace9bSHong Zhang     nzb = 0;    /* num of upper triangular entries in mat->B */
47216ebf90aSShri Abhyankar     for (i=0; i<m; i++) {
473e0bace9bSHong Zhang       nza   += (ai[i+1] - adiag[i]);
47416ebf90aSShri Abhyankar       countB = bi[i+1] - bi[i];
47516ebf90aSShri Abhyankar       bjj    = bj + bi[i];
476e0bace9bSHong Zhang       for (j=0; j<countB; j++) {
477e0bace9bSHong Zhang         if (garray[bjj[j]] > rstart) nzb++;
478e0bace9bSHong Zhang       }
479e0bace9bSHong Zhang     }
48016ebf90aSShri Abhyankar 
481e0bace9bSHong Zhang     nz   = nza + nzb; /* total nz of upper triangular part of mat */
48216ebf90aSShri Abhyankar     *nnz = nz;
483185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
484185f6596SHong Zhang     col  = row + nz;
485185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
486185f6596SHong Zhang 
48716ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
48816ebf90aSShri Abhyankar   } else {
48916ebf90aSShri Abhyankar     row = *r; col = *c; val = *v;
49016ebf90aSShri Abhyankar   }
49116ebf90aSShri Abhyankar 
49216ebf90aSShri Abhyankar   jj = 0; irow = rstart;
49316ebf90aSShri Abhyankar   for (i=0; i<m; i++) {
49416ebf90aSShri Abhyankar     ajj    = aj + adiag[i];                 /* ptr to the beginning of the diagonal of this row */
49516ebf90aSShri Abhyankar     v1     = av + adiag[i];
49616ebf90aSShri Abhyankar     countA = ai[i+1] - adiag[i];
49716ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
49816ebf90aSShri Abhyankar     bjj    = bj + bi[i];
49916ebf90aSShri Abhyankar     v2     = bv + bi[i];
50016ebf90aSShri Abhyankar 
50116ebf90aSShri Abhyankar     /* A-part */
50216ebf90aSShri Abhyankar     for (j=0; j<countA; j++) {
503bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
50416ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
50516ebf90aSShri Abhyankar       }
50616ebf90aSShri Abhyankar       val[jj++] = v1[j];
50716ebf90aSShri Abhyankar     }
50816ebf90aSShri Abhyankar 
50916ebf90aSShri Abhyankar     /* B-part */
51016ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
51116ebf90aSShri Abhyankar       if (garray[bjj[j]] > rstart) {
512bccb9932SShri Abhyankar         if (reuse == MAT_INITIAL_MATRIX) {
51316ebf90aSShri Abhyankar           row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
51416ebf90aSShri Abhyankar         }
51516ebf90aSShri Abhyankar         val[jj++] = v2[j];
51616ebf90aSShri Abhyankar       }
517397b6df1SKris Buschelman     }
518397b6df1SKris Buschelman     irow++;
519397b6df1SKris Buschelman   }
520397b6df1SKris Buschelman   PetscFunctionReturn(0);
521397b6df1SKris Buschelman }
522397b6df1SKris Buschelman 
523397b6df1SKris Buschelman #undef __FUNCT__
5243924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS"
525dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A)
526dfbe8321SBarry Smith {
527a5e57a09SHong Zhang   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
528dfbe8321SBarry Smith   PetscErrorCode ierr;
529b24902e0SBarry Smith 
530397b6df1SKris Buschelman   PetscFunctionBegin;
531a5e57a09SHong Zhang   if (mumps->CleanUpMUMPS) {
532397b6df1SKris Buschelman     /* Terminate instance, deallocate memories */
533a5e57a09SHong Zhang     ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr);
534a5e57a09SHong Zhang     ierr = VecScatterDestroy(&mumps->scat_rhs);CHKERRQ(ierr);
535a5e57a09SHong Zhang     ierr = VecDestroy(&mumps->b_seq);CHKERRQ(ierr);
536a5e57a09SHong Zhang     ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
537a5e57a09SHong Zhang     ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr);
538a5e57a09SHong Zhang     ierr = PetscFree(mumps->id.perm_in);CHKERRQ(ierr);
539a5e57a09SHong Zhang     ierr = PetscFree(mumps->irn);CHKERRQ(ierr);
5402205254eSKarl Rupp 
541a5e57a09SHong Zhang     mumps->id.job = JOB_END;
542a5e57a09SHong Zhang     PetscMUMPS_c(&mumps->id);
543a5e57a09SHong Zhang     ierr = MPI_Comm_free(&(mumps->comm_mumps));CHKERRQ(ierr);
544397b6df1SKris Buschelman   }
545a5e57a09SHong Zhang   if (mumps->Destroy) {
546a5e57a09SHong Zhang     ierr = (mumps->Destroy)(A);CHKERRQ(ierr);
547bf0cc555SLisandro Dalcin   }
548bf0cc555SLisandro Dalcin   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
549bf0cc555SLisandro Dalcin 
55097969023SHong Zhang   /* clear composed functions */
551bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr);
552bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetIcntl_C",NULL);CHKERRQ(ierr);
553bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetIcntl_C",NULL);CHKERRQ(ierr);
554bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetCntl_C",NULL);CHKERRQ(ierr);
555bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetCntl_C",NULL);CHKERRQ(ierr);
556bc6112feSHong Zhang 
557ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfo_C",NULL);CHKERRQ(ierr);
558ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfog_C",NULL);CHKERRQ(ierr);
559ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfo_C",NULL);CHKERRQ(ierr);
560ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfog_C",NULL);CHKERRQ(ierr);
561397b6df1SKris Buschelman   PetscFunctionReturn(0);
562397b6df1SKris Buschelman }
563397b6df1SKris Buschelman 
564397b6df1SKris Buschelman #undef __FUNCT__
565f6c57405SHong Zhang #define __FUNCT__ "MatSolve_MUMPS"
566b24902e0SBarry Smith PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x)
567b24902e0SBarry Smith {
568a5e57a09SHong Zhang   Mat_MUMPS        *mumps=(Mat_MUMPS*)A->spptr;
569d54de34fSKris Buschelman   PetscScalar      *array;
57067877ebaSShri Abhyankar   Vec              b_seq;
571329ec9b3SHong Zhang   IS               is_iden,is_petsc;
572dfbe8321SBarry Smith   PetscErrorCode   ierr;
573329ec9b3SHong Zhang   PetscInt         i;
574883f2eb9SBarry Smith   static PetscBool cite1 = PETSC_FALSE,cite2 = PETSC_FALSE;
575397b6df1SKris Buschelman 
576397b6df1SKris Buschelman   PetscFunctionBegin;
577883f2eb9SBarry 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);
578883f2eb9SBarry 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);
579a5e57a09SHong Zhang   mumps->id.nrhs = 1;
580a5e57a09SHong Zhang   b_seq          = mumps->b_seq;
581a5e57a09SHong Zhang   if (mumps->size > 1) {
582329ec9b3SHong Zhang     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
583a5e57a09SHong Zhang     ierr = VecScatterBegin(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
584a5e57a09SHong Zhang     ierr = VecScatterEnd(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
585a5e57a09SHong Zhang     if (!mumps->myid) {ierr = VecGetArray(b_seq,&array);CHKERRQ(ierr);}
586397b6df1SKris Buschelman   } else {  /* size == 1 */
587397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
588397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
589397b6df1SKris Buschelman   }
590a5e57a09SHong Zhang   if (!mumps->myid) { /* define rhs on the host */
591a5e57a09SHong Zhang     mumps->id.nrhs = 1;
592397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
5932907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
594a5e57a09SHong Zhang     mumps->id.rhs = (mumps_complex*)array;
5952907cef9SHong Zhang #else
596a5e57a09SHong Zhang     mumps->id.rhs = (mumps_double_complex*)array;
5972907cef9SHong Zhang #endif
598397b6df1SKris Buschelman #else
599a5e57a09SHong Zhang     mumps->id.rhs = array;
600397b6df1SKris Buschelman #endif
601397b6df1SKris Buschelman   }
602397b6df1SKris Buschelman 
603397b6df1SKris Buschelman   /* solve phase */
604329ec9b3SHong Zhang   /*-------------*/
605a5e57a09SHong Zhang   mumps->id.job = JOB_SOLVE;
606a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
607a5e57a09SHong 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));
608397b6df1SKris Buschelman 
609a5e57a09SHong Zhang   if (mumps->size > 1) { /* convert mumps distributed solution to petsc mpi x */
610a5e57a09SHong Zhang     if (mumps->scat_sol && mumps->ICNTL9_pre != mumps->id.ICNTL(9)) {
611a5e57a09SHong Zhang       /* when id.ICNTL(9) changes, the contents of lsol_loc may change (not its size, lsol_loc), recreates scat_sol */
612a5e57a09SHong Zhang       ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
613397b6df1SKris Buschelman     }
614a5e57a09SHong Zhang     if (!mumps->scat_sol) { /* create scatter scat_sol */
615a5e57a09SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,mumps->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */
616a5e57a09SHong Zhang       for (i=0; i<mumps->id.lsol_loc; i++) {
617a5e57a09SHong Zhang         mumps->id.isol_loc[i] -= 1; /* change Fortran style to C style */
618a5e57a09SHong Zhang       }
619a5e57a09SHong Zhang       ierr = ISCreateGeneral(PETSC_COMM_SELF,mumps->id.lsol_loc,mumps->id.isol_loc,PETSC_COPY_VALUES,&is_petsc);CHKERRQ(ierr);  /* to */
620a5e57a09SHong Zhang       ierr = VecScatterCreate(mumps->x_seq,is_iden,x,is_petsc,&mumps->scat_sol);CHKERRQ(ierr);
6216bf464f9SBarry Smith       ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
6226bf464f9SBarry Smith       ierr = ISDestroy(&is_petsc);CHKERRQ(ierr);
6232205254eSKarl Rupp 
624a5e57a09SHong Zhang       mumps->ICNTL9_pre = mumps->id.ICNTL(9); /* save current value of id.ICNTL(9) */
625397b6df1SKris Buschelman     }
626a5e57a09SHong Zhang 
627a5e57a09SHong Zhang     ierr = VecScatterBegin(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
628a5e57a09SHong Zhang     ierr = VecScatterEnd(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
629329ec9b3SHong Zhang   }
630397b6df1SKris Buschelman   PetscFunctionReturn(0);
631397b6df1SKris Buschelman }
632397b6df1SKris Buschelman 
63351d5961aSHong Zhang #undef __FUNCT__
63451d5961aSHong Zhang #define __FUNCT__ "MatSolveTranspose_MUMPS"
63551d5961aSHong Zhang PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x)
63651d5961aSHong Zhang {
637a5e57a09SHong Zhang   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
63851d5961aSHong Zhang   PetscErrorCode ierr;
63951d5961aSHong Zhang 
64051d5961aSHong Zhang   PetscFunctionBegin;
641a5e57a09SHong Zhang   mumps->id.ICNTL(9) = 0;
6422205254eSKarl Rupp 
6430ad0caddSJed Brown   ierr = MatSolve_MUMPS(A,b,x);CHKERRQ(ierr);
6442205254eSKarl Rupp 
645a5e57a09SHong Zhang   mumps->id.ICNTL(9) = 1;
64651d5961aSHong Zhang   PetscFunctionReturn(0);
64751d5961aSHong Zhang }
64851d5961aSHong Zhang 
649e0b74bf9SHong Zhang #undef __FUNCT__
650e0b74bf9SHong Zhang #define __FUNCT__ "MatMatSolve_MUMPS"
651e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X)
652e0b74bf9SHong Zhang {
653bda8bf91SBarry Smith   PetscErrorCode ierr;
654bda8bf91SBarry Smith   PetscBool      flg;
655bda8bf91SBarry Smith 
656e0b74bf9SHong Zhang   PetscFunctionBegin;
6570298fd71SBarry Smith   ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
658ce94432eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
6590298fd71SBarry Smith   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
660ce94432eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
6612205254eSKarl Rupp   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_MUMPS() is not implemented yet");
662e0b74bf9SHong Zhang   PetscFunctionReturn(0);
663e0b74bf9SHong Zhang }
664e0b74bf9SHong Zhang 
665ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX)
666a58c3f20SHong Zhang /*
667a58c3f20SHong Zhang   input:
668a58c3f20SHong Zhang    F:        numeric factor
669a58c3f20SHong Zhang   output:
670a58c3f20SHong Zhang    nneg:     total number of negative pivots
671a58c3f20SHong Zhang    nzero:    0
672a58c3f20SHong Zhang    npos:     (global dimension of F) - nneg
673a58c3f20SHong Zhang */
674a58c3f20SHong Zhang 
675a58c3f20SHong Zhang #undef __FUNCT__
676a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
677dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
678a58c3f20SHong Zhang {
679a5e57a09SHong Zhang   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
680dfbe8321SBarry Smith   PetscErrorCode ierr;
681c1490034SHong Zhang   PetscMPIInt    size;
682a58c3f20SHong Zhang 
683a58c3f20SHong Zhang   PetscFunctionBegin;
684ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)F),&size);CHKERRQ(ierr);
685bcb30aebSHong 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 */
686a5e57a09SHong 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));
687a58c3f20SHong Zhang   if (nneg) {
688a5e57a09SHong Zhang     if (!mumps->myid) {
689a5e57a09SHong Zhang       *nneg = mumps->id.INFOG(12);
690a58c3f20SHong Zhang     }
691a5e57a09SHong Zhang     ierr = MPI_Bcast(nneg,1,MPI_INT,0,mumps->comm_mumps);CHKERRQ(ierr);
692a58c3f20SHong Zhang   }
693a58c3f20SHong Zhang   if (nzero) *nzero = 0;
694d0f46423SBarry Smith   if (npos)  *npos  = F->rmap->N - (*nneg);
695a58c3f20SHong Zhang   PetscFunctionReturn(0);
696a58c3f20SHong Zhang }
697ace3df97SHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */
698a58c3f20SHong Zhang 
699397b6df1SKris Buschelman #undef __FUNCT__
700f6c57405SHong Zhang #define __FUNCT__ "MatFactorNumeric_MUMPS"
7010481f469SBarry Smith PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
702af281ebdSHong Zhang {
703a5e57a09SHong Zhang   Mat_MUMPS      *mumps =(Mat_MUMPS*)(F)->spptr;
7046849ba73SBarry Smith   PetscErrorCode ierr;
705e09efc27SHong Zhang   Mat            F_diag;
706ace3abfcSBarry Smith   PetscBool      isMPIAIJ;
707397b6df1SKris Buschelman 
708397b6df1SKris Buschelman   PetscFunctionBegin;
709a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
710397b6df1SKris Buschelman 
711397b6df1SKris Buschelman   /* numerical factorization phase */
712329ec9b3SHong Zhang   /*-------------------------------*/
713a5e57a09SHong Zhang   mumps->id.job = JOB_FACTNUMERIC;
714a5e57a09SHong Zhang   if (!mumps->id.ICNTL(18)) {
715a5e57a09SHong Zhang     if (!mumps->myid) {
716397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
7172907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
718a5e57a09SHong Zhang       mumps->id.a = (mumps_complex*)mumps->val;
7192907cef9SHong Zhang #else
720a5e57a09SHong Zhang       mumps->id.a = (mumps_double_complex*)mumps->val;
7212907cef9SHong Zhang #endif
722397b6df1SKris Buschelman #else
723a5e57a09SHong Zhang       mumps->id.a = mumps->val;
724397b6df1SKris Buschelman #endif
725397b6df1SKris Buschelman     }
726397b6df1SKris Buschelman   } else {
727397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
7282907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
729a5e57a09SHong Zhang     mumps->id.a_loc = (mumps_complex*)mumps->val;
7302907cef9SHong Zhang #else
731a5e57a09SHong Zhang     mumps->id.a_loc = (mumps_double_complex*)mumps->val;
7322907cef9SHong Zhang #endif
733397b6df1SKris Buschelman #else
734a5e57a09SHong Zhang     mumps->id.a_loc = mumps->val;
735397b6df1SKris Buschelman #endif
736397b6df1SKris Buschelman   }
737a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
738a5e57a09SHong Zhang   if (mumps->id.INFOG(1) < 0) {
739151787a6SHong Zhang     if (mumps->id.INFO(1) == -13) {
740151787a6SHong Zhang       if (mumps->id.INFO(2) < 0) {
741151787a6SHong 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));
742151787a6SHong Zhang       } else {
743151787a6SHong 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));
744151787a6SHong Zhang       }
745151787a6SHong 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));
746397b6df1SKris Buschelman   }
747a5e57a09SHong 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));
748397b6df1SKris Buschelman 
749a5e57a09SHong Zhang   if (mumps->size > 1) {
750251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr);
7512205254eSKarl Rupp     if (isMPIAIJ) F_diag = ((Mat_MPIAIJ*)(F)->data)->A;
7522205254eSKarl Rupp     else F_diag = ((Mat_MPISBAIJ*)(F)->data)->A;
753e09efc27SHong Zhang     F_diag->assembled = PETSC_TRUE;
754a5e57a09SHong Zhang     if (mumps->scat_sol) {
755a5e57a09SHong Zhang       ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
756a5e57a09SHong Zhang       ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr);
757a5e57a09SHong Zhang       ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr);
758329ec9b3SHong Zhang     }
7598ada1bb4SHong Zhang   }
760dcd589f8SShri Abhyankar   (F)->assembled      = PETSC_TRUE;
761a5e57a09SHong Zhang   mumps->matstruc     = SAME_NONZERO_PATTERN;
762a5e57a09SHong Zhang   mumps->CleanUpMUMPS = PETSC_TRUE;
76367877ebaSShri Abhyankar 
764a5e57a09SHong Zhang   if (mumps->size > 1) {
76567877ebaSShri Abhyankar     /* distributed solution */
766a5e57a09SHong Zhang     if (!mumps->scat_sol) {
76767877ebaSShri Abhyankar       /* Create x_seq=sol_loc for repeated use */
76867877ebaSShri Abhyankar       PetscInt    lsol_loc;
76967877ebaSShri Abhyankar       PetscScalar *sol_loc;
7702205254eSKarl Rupp 
771a5e57a09SHong Zhang       lsol_loc = mumps->id.INFO(23); /* length of sol_loc */
7722205254eSKarl Rupp 
773dcca6d9dSJed Brown       ierr = PetscMalloc2(lsol_loc,&sol_loc,lsol_loc,&mumps->id.isol_loc);CHKERRQ(ierr);
7742205254eSKarl Rupp 
775a5e57a09SHong Zhang       mumps->id.lsol_loc = lsol_loc;
77667877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
7772907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
778a5e57a09SHong Zhang       mumps->id.sol_loc = (mumps_complex*)sol_loc;
7792907cef9SHong Zhang #else
780a5e57a09SHong Zhang       mumps->id.sol_loc = (mumps_double_complex*)sol_loc;
7812907cef9SHong Zhang #endif
78267877ebaSShri Abhyankar #else
783a5e57a09SHong Zhang       mumps->id.sol_loc = sol_loc;
78467877ebaSShri Abhyankar #endif
785a5e57a09SHong Zhang       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq);CHKERRQ(ierr);
78667877ebaSShri Abhyankar     }
78767877ebaSShri Abhyankar   }
788397b6df1SKris Buschelman   PetscFunctionReturn(0);
789397b6df1SKris Buschelman }
790397b6df1SKris Buschelman 
7919a2535b5SHong Zhang /* Sets MUMPS options from the options database */
792dcd589f8SShri Abhyankar #undef __FUNCT__
7939a2535b5SHong Zhang #define __FUNCT__ "PetscSetMUMPSFromOptions"
7949a2535b5SHong Zhang PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A)
795dcd589f8SShri Abhyankar {
7969a2535b5SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
797dcd589f8SShri Abhyankar   PetscErrorCode ierr;
798dcd589f8SShri Abhyankar   PetscInt       icntl;
799ace3abfcSBarry Smith   PetscBool      flg;
800dcd589f8SShri Abhyankar 
801dcd589f8SShri Abhyankar   PetscFunctionBegin;
802ce94432eSBarry Smith   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
8039a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);CHKERRQ(ierr);
8049a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(1) = icntl;
8059a2535b5SHong 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);
8069a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(2) = icntl;
8079a2535b5SHong 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);
8089a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(3) = icntl;
809dcd589f8SShri Abhyankar 
8109a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",mumps->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
8119a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(4) = icntl;
8129a2535b5SHong Zhang   if (mumps->id.ICNTL(4) || PetscLogPrintInfo) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */
8139a2535b5SHong Zhang 
8149a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_6","ICNTL(6): permuting and/or scaling the matrix (0 to 7)","None",mumps->id.ICNTL(6),&icntl,&flg);CHKERRQ(ierr);
8159a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(6) = icntl;
8169a2535b5SHong Zhang 
8179a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7). 3=Scotch, 4=PORD, 5=Metis","None",mumps->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
818dcd589f8SShri Abhyankar   if (flg) {
8192205254eSKarl 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");
8202205254eSKarl Rupp     else mumps->id.ICNTL(7) = icntl;
821dcd589f8SShri Abhyankar   }
822e0b74bf9SHong Zhang 
8230298fd71SBarry 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);
8240298fd71SBarry 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);
8250298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_11","ICNTL(11): statistics related to the linear system solved (via -ksp_view)","None",mumps->id.ICNTL(11),&mumps->id.ICNTL(11),NULL);CHKERRQ(ierr);
8260298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control: defines the ordering strategy with scaling constraints (0 to 3)","None",mumps->id.ICNTL(12),&mumps->id.ICNTL(12),NULL);CHKERRQ(ierr);
8270298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control: with or without ScaLAPACK","None",mumps->id.ICNTL(13),&mumps->id.ICNTL(13),NULL);CHKERRQ(ierr);
8280298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage of estimated workspace increase","None",mumps->id.ICNTL(14),&mumps->id.ICNTL(14),NULL);CHKERRQ(ierr);
8290298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): Schur complement","None",mumps->id.ICNTL(19),&mumps->id.ICNTL(19),NULL);CHKERRQ(ierr);
8309a2535b5SHong Zhang 
8310298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_22","ICNTL(22): in-core/out-of-core facility (0 or 1)","None",mumps->id.ICNTL(22),&mumps->id.ICNTL(22),NULL);CHKERRQ(ierr);
8320298fd71SBarry 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);
8330298fd71SBarry 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);
8349a2535b5SHong Zhang   if (mumps->id.ICNTL(24)) {
8359a2535b5SHong Zhang     mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */
836d7ebd59bSHong Zhang   }
837d7ebd59bSHong Zhang 
8380298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_25","ICNTL(25): computation of a null space basis","None",mumps->id.ICNTL(25),&mumps->id.ICNTL(25),NULL);CHKERRQ(ierr);
8390298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_26","ICNTL(26): Schur options for right-hand side or solution vector","None",mumps->id.ICNTL(26),&mumps->id.ICNTL(26),NULL);CHKERRQ(ierr);
8400298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_27","ICNTL(27): experimental parameter","None",mumps->id.ICNTL(27),&mumps->id.ICNTL(27),NULL);CHKERRQ(ierr);
8410298fd71SBarry 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);
8420298fd71SBarry Smith   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);
8430298fd71SBarry 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);
8440298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_31","ICNTL(31): factors can be discarded in the solve phase","None",mumps->id.ICNTL(31),&mumps->id.ICNTL(31),NULL);CHKERRQ(ierr);
8450298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),NULL);CHKERRQ(ierr);
846dcd589f8SShri Abhyankar 
8470298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),NULL);CHKERRQ(ierr);
8480298fd71SBarry 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);
8490298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),NULL);CHKERRQ(ierr);
8500298fd71SBarry 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);
8510298fd71SBarry 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);
852e5bb22a1SHong Zhang 
8530298fd71SBarry Smith   ierr = PetscOptionsString("-mat_mumps_ooc_tmpdir", "out of core directory", "None", mumps->id.ooc_tmpdir, mumps->id.ooc_tmpdir, 256, NULL);
854dcd589f8SShri Abhyankar   PetscOptionsEnd();
855dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
856dcd589f8SShri Abhyankar }
857dcd589f8SShri Abhyankar 
858dcd589f8SShri Abhyankar #undef __FUNCT__
859dcd589f8SShri Abhyankar #define __FUNCT__ "PetscInitializeMUMPS"
860f697e70eSHong Zhang PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps)
861dcd589f8SShri Abhyankar {
862dcd589f8SShri Abhyankar   PetscErrorCode ierr;
863dcd589f8SShri Abhyankar 
864dcd589f8SShri Abhyankar   PetscFunctionBegin;
865ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)A), &mumps->myid);
866ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&mumps->size);CHKERRQ(ierr);
867ce94432eSBarry Smith   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(mumps->comm_mumps));CHKERRQ(ierr);
8682205254eSKarl Rupp 
869f697e70eSHong Zhang   mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps);
870f697e70eSHong Zhang 
871f697e70eSHong Zhang   mumps->id.job = JOB_INIT;
872f697e70eSHong Zhang   mumps->id.par = 1;  /* host participates factorizaton and solve */
873f697e70eSHong Zhang   mumps->id.sym = mumps->sym;
8742907cef9SHong Zhang   PetscMUMPS_c(&mumps->id);
875f697e70eSHong Zhang 
876f697e70eSHong Zhang   mumps->CleanUpMUMPS = PETSC_FALSE;
8770298fd71SBarry Smith   mumps->scat_rhs     = NULL;
8780298fd71SBarry Smith   mumps->scat_sol     = NULL;
8799a2535b5SHong Zhang 
88070544d5fSHong Zhang   /* set PETSc-MUMPS default options - override MUMPS default */
8819a2535b5SHong Zhang   mumps->id.ICNTL(3) = 0;
8829a2535b5SHong Zhang   mumps->id.ICNTL(4) = 0;
8839a2535b5SHong Zhang   if (mumps->size == 1) {
8849a2535b5SHong Zhang     mumps->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
8859a2535b5SHong Zhang   } else {
8869a2535b5SHong Zhang     mumps->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
88770544d5fSHong Zhang     mumps->id.ICNTL(21) = 1;   /* distributed solution */
8889a2535b5SHong Zhang   }
889dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
890dcd589f8SShri Abhyankar }
891dcd589f8SShri Abhyankar 
892a5e57a09SHong 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 */
893397b6df1SKris Buschelman #undef __FUNCT__
894f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
8950481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
896b24902e0SBarry Smith {
897a5e57a09SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
898dcd589f8SShri Abhyankar   PetscErrorCode ierr;
89967877ebaSShri Abhyankar   Vec            b;
90067877ebaSShri Abhyankar   IS             is_iden;
90167877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
902397b6df1SKris Buschelman 
903397b6df1SKris Buschelman   PetscFunctionBegin;
904a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
905dcd589f8SShri Abhyankar 
9069a2535b5SHong Zhang   /* Set MUMPS options from the options database */
9079a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
908dcd589f8SShri Abhyankar 
909a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
910dcd589f8SShri Abhyankar 
91167877ebaSShri Abhyankar   /* analysis phase */
91267877ebaSShri Abhyankar   /*----------------*/
913a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
914a5e57a09SHong Zhang   mumps->id.n   = M;
915a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
91667877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
917a5e57a09SHong Zhang     if (!mumps->myid) {
918a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
919a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
92067877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
9212907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
922a5e57a09SHong Zhang         mumps->id.a = (mumps_complex*)mumps->val;
9232907cef9SHong Zhang #else
924a5e57a09SHong Zhang         mumps->id.a = (mumps_double_complex*)mumps->val;
9252907cef9SHong Zhang #endif
92667877ebaSShri Abhyankar #else
927a5e57a09SHong Zhang         mumps->id.a = mumps->val;
92867877ebaSShri Abhyankar #endif
92967877ebaSShri Abhyankar       }
930a5e57a09SHong Zhang       if (mumps->id.ICNTL(7) == 1) { /* use user-provide matrix ordering - assuming r = c ordering */
9315248a706SHong Zhang         /*
9325248a706SHong Zhang         PetscBool      flag;
9335248a706SHong Zhang         ierr = ISEqual(r,c,&flag);CHKERRQ(ierr);
9345248a706SHong Zhang         if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"row_perm != col_perm");
9355248a706SHong Zhang         ierr = ISView(r,PETSC_VIEWER_STDOUT_SELF);
9365248a706SHong Zhang          */
937a5e57a09SHong Zhang         if (!mumps->myid) {
938e0b74bf9SHong Zhang           const PetscInt *idx;
939e0b74bf9SHong Zhang           PetscInt       i,*perm_in;
9402205254eSKarl Rupp 
941785e854fSJed Brown           ierr = PetscMalloc1(M,&perm_in);CHKERRQ(ierr);
942e0b74bf9SHong Zhang           ierr = ISGetIndices(r,&idx);CHKERRQ(ierr);
9432205254eSKarl Rupp 
944a5e57a09SHong Zhang           mumps->id.perm_in = perm_in;
945e0b74bf9SHong Zhang           for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */
946e0b74bf9SHong Zhang           ierr = ISRestoreIndices(r,&idx);CHKERRQ(ierr);
947e0b74bf9SHong Zhang         }
948e0b74bf9SHong Zhang       }
94967877ebaSShri Abhyankar     }
95067877ebaSShri Abhyankar     break;
95167877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
952a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
953a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
954a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
95567877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
9562907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
957a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_complex*)mumps->val;
9582907cef9SHong Zhang #else
959a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_double_complex*)mumps->val;
9602907cef9SHong Zhang #endif
96167877ebaSShri Abhyankar #else
962a5e57a09SHong Zhang       mumps->id.a_loc = mumps->val;
96367877ebaSShri Abhyankar #endif
96467877ebaSShri Abhyankar     }
96567877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
966a5e57a09SHong Zhang     if (!mumps->myid) {
967a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
96867877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
96967877ebaSShri Abhyankar     } else {
970a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
97167877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
97267877ebaSShri Abhyankar     }
973c0dedaeaSBarry Smith     ierr = MatGetVecs(A,NULL,&b);CHKERRQ(ierr);
974a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
9756bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
9766bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
97767877ebaSShri Abhyankar     break;
97867877ebaSShri Abhyankar   }
979a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
980a5e57a09SHong 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));
98167877ebaSShri Abhyankar 
982719d5645SBarry Smith   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
983dcd589f8SShri Abhyankar   F->ops->solve           = MatSolve_MUMPS;
98451d5961aSHong Zhang   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
98517f96c7aSHong Zhang   F->ops->matsolve        = 0;  /* use MatMatSolve_Basic() until mumps supports distributed rhs */
986b24902e0SBarry Smith   PetscFunctionReturn(0);
987b24902e0SBarry Smith }
988b24902e0SBarry Smith 
989450b117fSShri Abhyankar /* Note the Petsc r and c permutations are ignored */
990450b117fSShri Abhyankar #undef __FUNCT__
991450b117fSShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_BAIJMUMPS"
992450b117fSShri Abhyankar PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
993450b117fSShri Abhyankar {
994a5e57a09SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
995dcd589f8SShri Abhyankar   PetscErrorCode ierr;
99667877ebaSShri Abhyankar   Vec            b;
99767877ebaSShri Abhyankar   IS             is_iden;
99867877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
999450b117fSShri Abhyankar 
1000450b117fSShri Abhyankar   PetscFunctionBegin;
1001a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1002dcd589f8SShri Abhyankar 
10039a2535b5SHong Zhang   /* Set MUMPS options from the options database */
10049a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1005dcd589f8SShri Abhyankar 
1006a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
100767877ebaSShri Abhyankar 
100867877ebaSShri Abhyankar   /* analysis phase */
100967877ebaSShri Abhyankar   /*----------------*/
1010a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1011a5e57a09SHong Zhang   mumps->id.n   = M;
1012a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
101367877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1014a5e57a09SHong Zhang     if (!mumps->myid) {
1015a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1016a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
101767877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
10182907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1019a5e57a09SHong Zhang         mumps->id.a = (mumps_complex*)mumps->val;
10202907cef9SHong Zhang #else
1021a5e57a09SHong Zhang         mumps->id.a = (mumps_double_complex*)mumps->val;
10222907cef9SHong Zhang #endif
102367877ebaSShri Abhyankar #else
1024a5e57a09SHong Zhang         mumps->id.a = mumps->val;
102567877ebaSShri Abhyankar #endif
102667877ebaSShri Abhyankar       }
102767877ebaSShri Abhyankar     }
102867877ebaSShri Abhyankar     break;
102967877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1030a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1031a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1032a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
103367877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
10342907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1035a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_complex*)mumps->val;
10362907cef9SHong Zhang #else
1037a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_double_complex*)mumps->val;
10382907cef9SHong Zhang #endif
103967877ebaSShri Abhyankar #else
1040a5e57a09SHong Zhang       mumps->id.a_loc = mumps->val;
104167877ebaSShri Abhyankar #endif
104267877ebaSShri Abhyankar     }
104367877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1044a5e57a09SHong Zhang     if (!mumps->myid) {
1045a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
104667877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
104767877ebaSShri Abhyankar     } else {
1048a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
104967877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
105067877ebaSShri Abhyankar     }
1051c0dedaeaSBarry Smith     ierr = MatGetVecs(A,NULL,&b);CHKERRQ(ierr);
1052a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
10536bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
10546bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
105567877ebaSShri Abhyankar     break;
105667877ebaSShri Abhyankar   }
1057a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
1058a5e57a09SHong 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));
105967877ebaSShri Abhyankar 
1060450b117fSShri Abhyankar   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1061dcd589f8SShri Abhyankar   F->ops->solve           = MatSolve_MUMPS;
106251d5961aSHong Zhang   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1063450b117fSShri Abhyankar   PetscFunctionReturn(0);
1064450b117fSShri Abhyankar }
1065b24902e0SBarry Smith 
1066141f4205SHong Zhang /* Note the Petsc r permutation and factor info are ignored */
1067397b6df1SKris Buschelman #undef __FUNCT__
106867877ebaSShri Abhyankar #define __FUNCT__ "MatCholeskyFactorSymbolic_MUMPS"
106967877ebaSShri Abhyankar PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
1070b24902e0SBarry Smith {
1071a5e57a09SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
1072dcd589f8SShri Abhyankar   PetscErrorCode ierr;
107367877ebaSShri Abhyankar   Vec            b;
107467877ebaSShri Abhyankar   IS             is_iden;
107567877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1076397b6df1SKris Buschelman 
1077397b6df1SKris Buschelman   PetscFunctionBegin;
1078a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1079dcd589f8SShri Abhyankar 
10809a2535b5SHong Zhang   /* Set MUMPS options from the options database */
10819a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1082dcd589f8SShri Abhyankar 
1083a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1084dcd589f8SShri Abhyankar 
108567877ebaSShri Abhyankar   /* analysis phase */
108667877ebaSShri Abhyankar   /*----------------*/
1087a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1088a5e57a09SHong Zhang   mumps->id.n   = M;
1089a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
109067877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1091a5e57a09SHong Zhang     if (!mumps->myid) {
1092a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1093a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
109467877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
10952907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1096a5e57a09SHong Zhang         mumps->id.a = (mumps_complex*)mumps->val;
10972907cef9SHong Zhang #else
1098a5e57a09SHong Zhang         mumps->id.a = (mumps_double_complex*)mumps->val;
10992907cef9SHong Zhang #endif
110067877ebaSShri Abhyankar #else
1101a5e57a09SHong Zhang         mumps->id.a = mumps->val;
110267877ebaSShri Abhyankar #endif
110367877ebaSShri Abhyankar       }
110467877ebaSShri Abhyankar     }
110567877ebaSShri Abhyankar     break;
110667877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1107a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1108a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1109a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
111067877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
11112907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1112a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_complex*)mumps->val;
11132907cef9SHong Zhang #else
1114a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_double_complex*)mumps->val;
11152907cef9SHong Zhang #endif
111667877ebaSShri Abhyankar #else
1117a5e57a09SHong Zhang       mumps->id.a_loc = mumps->val;
111867877ebaSShri Abhyankar #endif
111967877ebaSShri Abhyankar     }
112067877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1121a5e57a09SHong Zhang     if (!mumps->myid) {
1122a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
112367877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
112467877ebaSShri Abhyankar     } else {
1125a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
112667877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
112767877ebaSShri Abhyankar     }
1128c0dedaeaSBarry Smith     ierr = MatGetVecs(A,NULL,&b);CHKERRQ(ierr);
1129a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
11306bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
11316bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
113267877ebaSShri Abhyankar     break;
113367877ebaSShri Abhyankar   }
1134a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
1135a5e57a09SHong 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));
113667877ebaSShri Abhyankar 
11372792810eSHong Zhang   F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
1138dcd589f8SShri Abhyankar   F->ops->solve                 = MatSolve_MUMPS;
113951d5961aSHong Zhang   F->ops->solvetranspose        = MatSolve_MUMPS;
114030c107b7SHong Zhang   F->ops->matsolve              = 0; /* use MatMatSolve_Basic() until mumps supports distributed rhs */
1141db4efbfdSBarry Smith #if !defined(PETSC_USE_COMPLEX)
114205aa0992SJose Roman   F->ops->getinertia = MatGetInertia_SBAIJMUMPS;
114305aa0992SJose Roman #else
11440298fd71SBarry Smith   F->ops->getinertia = NULL;
1145db4efbfdSBarry Smith #endif
1146b24902e0SBarry Smith   PetscFunctionReturn(0);
1147b24902e0SBarry Smith }
1148b24902e0SBarry Smith 
1149397b6df1SKris Buschelman #undef __FUNCT__
115064e6c443SBarry Smith #define __FUNCT__ "MatView_MUMPS"
115164e6c443SBarry Smith PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
115274ed9c26SBarry Smith {
1153f6c57405SHong Zhang   PetscErrorCode    ierr;
115464e6c443SBarry Smith   PetscBool         iascii;
115564e6c443SBarry Smith   PetscViewerFormat format;
1156a5e57a09SHong Zhang   Mat_MUMPS         *mumps=(Mat_MUMPS*)A->spptr;
1157f6c57405SHong Zhang 
1158f6c57405SHong Zhang   PetscFunctionBegin;
115964e6c443SBarry Smith   /* check if matrix is mumps type */
116064e6c443SBarry Smith   if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0);
116164e6c443SBarry Smith 
1162251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
116364e6c443SBarry Smith   if (iascii) {
116464e6c443SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
116564e6c443SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO) {
116664e6c443SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
1167a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                   %d \n",mumps->id.sym);CHKERRQ(ierr);
1168a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):            %d \n",mumps->id.par);CHKERRQ(ierr);
1169a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):         %d \n",mumps->id.ICNTL(1));CHKERRQ(ierr);
1170a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg): %d \n",mumps->id.ICNTL(2));CHKERRQ(ierr);
1171a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):   %d \n",mumps->id.ICNTL(3));CHKERRQ(ierr);
1172a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):        %d \n",mumps->id.ICNTL(4));CHKERRQ(ierr);
1173a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):         %d \n",mumps->id.ICNTL(5));CHKERRQ(ierr);
1174a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):        %d \n",mumps->id.ICNTL(6));CHKERRQ(ierr);
1175a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (sequentia matrix ordering):%d \n",mumps->id.ICNTL(7));CHKERRQ(ierr);
1176a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scalling strategy):        %d \n",mumps->id.ICNTL(8));CHKERRQ(ierr);
1177a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements):  %d \n",mumps->id.ICNTL(10));CHKERRQ(ierr);
1178a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):          %d \n",mumps->id.ICNTL(11));CHKERRQ(ierr);
1179a5e57a09SHong Zhang       if (mumps->id.ICNTL(11)>0) {
1180a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(4) (inf norm of input mat):        %g\n",mumps->id.RINFOG(4));CHKERRQ(ierr);
1181a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(5) (inf norm of solution):         %g\n",mumps->id.RINFOG(5));CHKERRQ(ierr);
1182a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(6) (inf norm of residual):         %g\n",mumps->id.RINFOG(6));CHKERRQ(ierr);
1183a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8));CHKERRQ(ierr);
1184a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(9) (error estimate):               %g \n",mumps->id.RINFOG(9));CHKERRQ(ierr);
1185a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11));CHKERRQ(ierr);
1186f6c57405SHong Zhang       }
1187a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",mumps->id.ICNTL(12));CHKERRQ(ierr);
1188a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",mumps->id.ICNTL(13));CHKERRQ(ierr);
1189a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",mumps->id.ICNTL(14));CHKERRQ(ierr);
1190f6c57405SHong Zhang       /* ICNTL(15-17) not used */
1191a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",mumps->id.ICNTL(18));CHKERRQ(ierr);
1192a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Shur complement info):                       %d \n",mumps->id.ICNTL(19));CHKERRQ(ierr);
1193a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",mumps->id.ICNTL(20));CHKERRQ(ierr);
1194a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (somumpstion struct):                            %d \n",mumps->id.ICNTL(21));CHKERRQ(ierr);
1195a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",mumps->id.ICNTL(22));CHKERRQ(ierr);
1196a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",mumps->id.ICNTL(23));CHKERRQ(ierr);
1197c0165424SHong Zhang 
1198a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",mumps->id.ICNTL(24));CHKERRQ(ierr);
1199a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",mumps->id.ICNTL(25));CHKERRQ(ierr);
1200a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",mumps->id.ICNTL(26));CHKERRQ(ierr);
1201a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",mumps->id.ICNTL(27));CHKERRQ(ierr);
1202a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(28) (use parallel or sequential ordering):        %d \n",mumps->id.ICNTL(28));CHKERRQ(ierr);
1203a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(29) (parallel ordering):                          %d \n",mumps->id.ICNTL(29));CHKERRQ(ierr);
120442179a6aSHong Zhang 
1205a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(30) (user-specified set of entries in inv(A)):    %d \n",mumps->id.ICNTL(30));CHKERRQ(ierr);
1206a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(31) (factors is discarded in the solve phase):    %d \n",mumps->id.ICNTL(31));CHKERRQ(ierr);
1207a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(33) (compute determinant):                        %d \n",mumps->id.ICNTL(33));CHKERRQ(ierr);
1208f6c57405SHong Zhang 
1209a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",mumps->id.CNTL(1));CHKERRQ(ierr);
1210a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2));CHKERRQ(ierr);
1211a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absomumpste pivoting threshold):      %g \n",mumps->id.CNTL(3));CHKERRQ(ierr);
1212a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(4) (vamumpse of static pivoting):         %g \n",mumps->id.CNTL(4));CHKERRQ(ierr);
1213a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",mumps->id.CNTL(5));CHKERRQ(ierr);
1214f6c57405SHong Zhang 
1215f6c57405SHong Zhang       /* infomation local to each processor */
121634ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);
12177b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
1218a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %g \n",mumps->myid,mumps->id.RINFO(1));CHKERRQ(ierr);
121934ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
122034ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);
1221a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(2));CHKERRQ(ierr);
122234ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
122334ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);
1224a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(3));CHKERRQ(ierr);
122534ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1226f6c57405SHong Zhang 
122734ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr);
1228a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] %d \n",mumps->myid,mumps->id.INFO(15));CHKERRQ(ierr);
122934ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1230f6c57405SHong Zhang 
123134ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr);
1232a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(16));CHKERRQ(ierr);
123334ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1234f6c57405SHong Zhang 
123534ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);
1236a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(23));CHKERRQ(ierr);
123734ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
12387b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
1239f6c57405SHong Zhang 
1240a5e57a09SHong Zhang       if (!mumps->myid) { /* information from the host */
1241a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",mumps->id.RINFOG(1));CHKERRQ(ierr);
1242a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",mumps->id.RINFOG(2));CHKERRQ(ierr);
1243a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",mumps->id.RINFOG(3));CHKERRQ(ierr);
1244a5e57a09SHong 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);
1245f6c57405SHong Zhang 
1246a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(3));CHKERRQ(ierr);
1247a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(4));CHKERRQ(ierr);
1248a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",mumps->id.INFOG(5));CHKERRQ(ierr);
1249a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",mumps->id.INFOG(6));CHKERRQ(ierr);
1250a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively use after analysis): %d \n",mumps->id.INFOG(7));CHKERRQ(ierr);
1251a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",mumps->id.INFOG(8));CHKERRQ(ierr);
1252a5e57a09SHong 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);
1253a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",mumps->id.INFOG(10));CHKERRQ(ierr);
1254a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",mumps->id.INFOG(11));CHKERRQ(ierr);
1255a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",mumps->id.INFOG(12));CHKERRQ(ierr);
1256a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",mumps->id.INFOG(13));CHKERRQ(ierr);
1257a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",mumps->id.INFOG(14));CHKERRQ(ierr);
1258a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",mumps->id.INFOG(15));CHKERRQ(ierr);
1259a5e57a09SHong 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);
1260a5e57a09SHong 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);
1261a5e57a09SHong 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);
1262a5e57a09SHong 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);
1263a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",mumps->id.INFOG(20));CHKERRQ(ierr);
1264a5e57a09SHong 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);
1265a5e57a09SHong 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);
1266a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",mumps->id.INFOG(23));CHKERRQ(ierr);
1267a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",mumps->id.INFOG(24));CHKERRQ(ierr);
1268a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",mumps->id.INFOG(25));CHKERRQ(ierr);
126940d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(28) (after factorization: number of null pivots encountered): %d\n",mumps->id.INFOG(28));CHKERRQ(ierr);
127040d435e3SHong 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);
127140d435e3SHong 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);
127240d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(32) (after analysis: type of analysis done): %d\n",mumps->id.INFOG(32));CHKERRQ(ierr);
127340d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(33) (value used for ICNTL(8)): %d\n",mumps->id.INFOG(33));CHKERRQ(ierr);
127440d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(34) (exponent of the determinant if determinant is requested): %d\n",mumps->id.INFOG(34));CHKERRQ(ierr);
1275f6c57405SHong Zhang       }
1276f6c57405SHong Zhang     }
1277cb828f0fSHong Zhang   }
1278f6c57405SHong Zhang   PetscFunctionReturn(0);
1279f6c57405SHong Zhang }
1280f6c57405SHong Zhang 
128135bd34faSBarry Smith #undef __FUNCT__
128235bd34faSBarry Smith #define __FUNCT__ "MatGetInfo_MUMPS"
128335bd34faSBarry Smith PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
128435bd34faSBarry Smith {
1285cb828f0fSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)A->spptr;
128635bd34faSBarry Smith 
128735bd34faSBarry Smith   PetscFunctionBegin;
128835bd34faSBarry Smith   info->block_size        = 1.0;
1289cb828f0fSHong Zhang   info->nz_allocated      = mumps->id.INFOG(20);
1290cb828f0fSHong Zhang   info->nz_used           = mumps->id.INFOG(20);
129135bd34faSBarry Smith   info->nz_unneeded       = 0.0;
129235bd34faSBarry Smith   info->assemblies        = 0.0;
129335bd34faSBarry Smith   info->mallocs           = 0.0;
129435bd34faSBarry Smith   info->memory            = 0.0;
129535bd34faSBarry Smith   info->fill_ratio_given  = 0;
129635bd34faSBarry Smith   info->fill_ratio_needed = 0;
129735bd34faSBarry Smith   info->factor_mallocs    = 0;
129835bd34faSBarry Smith   PetscFunctionReturn(0);
129935bd34faSBarry Smith }
130035bd34faSBarry Smith 
13015ccb76cbSHong Zhang /* -------------------------------------------------------------------------------------------*/
13025ccb76cbSHong Zhang #undef __FUNCT__
13035ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl_MUMPS"
13045ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival)
13055ccb76cbSHong Zhang {
1306a5e57a09SHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
13075ccb76cbSHong Zhang 
13085ccb76cbSHong Zhang   PetscFunctionBegin;
1309a5e57a09SHong Zhang   mumps->id.ICNTL(icntl) = ival;
13105ccb76cbSHong Zhang   PetscFunctionReturn(0);
13115ccb76cbSHong Zhang }
13125ccb76cbSHong Zhang 
13135ccb76cbSHong Zhang #undef __FUNCT__
1314bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetIcntl_MUMPS"
1315bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt *ival)
1316bc6112feSHong Zhang {
1317bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1318bc6112feSHong Zhang 
1319bc6112feSHong Zhang   PetscFunctionBegin;
1320bc6112feSHong Zhang   *ival = mumps->id.ICNTL(icntl);
1321bc6112feSHong Zhang   PetscFunctionReturn(0);
1322bc6112feSHong Zhang }
1323bc6112feSHong Zhang 
1324bc6112feSHong Zhang #undef __FUNCT__
13255ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl"
13265ccb76cbSHong Zhang /*@
13275ccb76cbSHong Zhang   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()
13285ccb76cbSHong Zhang 
13295ccb76cbSHong Zhang    Logically Collective on Mat
13305ccb76cbSHong Zhang 
13315ccb76cbSHong Zhang    Input Parameters:
13325ccb76cbSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
13335ccb76cbSHong Zhang .  icntl - index of MUMPS parameter array ICNTL()
13345ccb76cbSHong Zhang -  ival - value of MUMPS ICNTL(icntl)
13355ccb76cbSHong Zhang 
13365ccb76cbSHong Zhang   Options Database:
13375ccb76cbSHong Zhang .   -mat_mumps_icntl_<icntl> <ival>
13385ccb76cbSHong Zhang 
13395ccb76cbSHong Zhang    Level: beginner
13405ccb76cbSHong Zhang 
13415ccb76cbSHong Zhang    References: MUMPS Users' Guide
13425ccb76cbSHong Zhang 
13435ccb76cbSHong Zhang .seealso: MatGetFactor()
13445ccb76cbSHong Zhang @*/
13455ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival)
13465ccb76cbSHong Zhang {
13475ccb76cbSHong Zhang   PetscErrorCode ierr;
13485ccb76cbSHong Zhang 
13495ccb76cbSHong Zhang   PetscFunctionBegin;
13505ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
13515ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,ival,3);
13525ccb76cbSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr);
13535ccb76cbSHong Zhang   PetscFunctionReturn(0);
13545ccb76cbSHong Zhang }
13555ccb76cbSHong Zhang 
1356bc6112feSHong Zhang #undef __FUNCT__
1357bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetIcntl"
1358*a21f80fcSHong Zhang /*@
1359*a21f80fcSHong Zhang   MatMumpsGetIcntl - Get MUMPS parameter ICNTL()
1360*a21f80fcSHong Zhang 
1361*a21f80fcSHong Zhang    Logically Collective on Mat
1362*a21f80fcSHong Zhang 
1363*a21f80fcSHong Zhang    Input Parameters:
1364*a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1365*a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array ICNTL()
1366*a21f80fcSHong Zhang 
1367*a21f80fcSHong Zhang   Output Parameter:
1368*a21f80fcSHong Zhang .  ival - value of MUMPS ICNTL(icntl)
1369*a21f80fcSHong Zhang 
1370*a21f80fcSHong Zhang    Level: beginner
1371*a21f80fcSHong Zhang 
1372*a21f80fcSHong Zhang    References: MUMPS Users' Guide
1373*a21f80fcSHong Zhang 
1374*a21f80fcSHong Zhang .seealso: MatGetFactor()
1375*a21f80fcSHong Zhang @*/
1376bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl(Mat F,PetscInt icntl,PetscInt *ival)
1377bc6112feSHong Zhang {
1378bc6112feSHong Zhang   PetscErrorCode ierr;
1379bc6112feSHong Zhang 
1380bc6112feSHong Zhang   PetscFunctionBegin;
1381bc6112feSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
1382bc6112feSHong Zhang   PetscValidIntPointer(ival,3);
1383bc6112feSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetIcntl_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
1384bc6112feSHong Zhang   PetscFunctionReturn(0);
1385bc6112feSHong Zhang }
1386bc6112feSHong Zhang 
13878928b65cSHong Zhang /* -------------------------------------------------------------------------------------------*/
13888928b65cSHong Zhang #undef __FUNCT__
13898928b65cSHong Zhang #define __FUNCT__ "MatMumpsSetCntl_MUMPS"
13908928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal val)
13918928b65cSHong Zhang {
13928928b65cSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
13938928b65cSHong Zhang 
13948928b65cSHong Zhang   PetscFunctionBegin;
13958928b65cSHong Zhang   mumps->id.CNTL(icntl) = val;
13968928b65cSHong Zhang   PetscFunctionReturn(0);
13978928b65cSHong Zhang }
13988928b65cSHong Zhang 
13998928b65cSHong Zhang #undef __FUNCT__
1400bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetCntl_MUMPS"
1401bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal *val)
1402bc6112feSHong Zhang {
1403bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1404bc6112feSHong Zhang 
1405bc6112feSHong Zhang   PetscFunctionBegin;
1406bc6112feSHong Zhang   *val = mumps->id.CNTL(icntl);
1407bc6112feSHong Zhang   PetscFunctionReturn(0);
1408bc6112feSHong Zhang }
1409bc6112feSHong Zhang 
1410bc6112feSHong Zhang #undef __FUNCT__
14118928b65cSHong Zhang #define __FUNCT__ "MatMumpsSetCntl"
14128928b65cSHong Zhang /*@
14138928b65cSHong Zhang   MatMumpsSetCntl - Set MUMPS parameter CNTL()
14148928b65cSHong Zhang 
14158928b65cSHong Zhang    Logically Collective on Mat
14168928b65cSHong Zhang 
14178928b65cSHong Zhang    Input Parameters:
14188928b65cSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
14198928b65cSHong Zhang .  icntl - index of MUMPS parameter array CNTL()
14208928b65cSHong Zhang -  val - value of MUMPS CNTL(icntl)
14218928b65cSHong Zhang 
14228928b65cSHong Zhang   Options Database:
14238928b65cSHong Zhang .   -mat_mumps_cntl_<icntl> <val>
14248928b65cSHong Zhang 
14258928b65cSHong Zhang    Level: beginner
14268928b65cSHong Zhang 
14278928b65cSHong Zhang    References: MUMPS Users' Guide
14288928b65cSHong Zhang 
14298928b65cSHong Zhang .seealso: MatGetFactor()
14308928b65cSHong Zhang @*/
14318928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl(Mat F,PetscInt icntl,PetscReal val)
14328928b65cSHong Zhang {
14338928b65cSHong Zhang   PetscErrorCode ierr;
14348928b65cSHong Zhang 
14358928b65cSHong Zhang   PetscFunctionBegin;
14368928b65cSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
1437bc6112feSHong Zhang   PetscValidLogicalCollectiveReal(F,val,3);
14388928b65cSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsSetCntl_C",(Mat,PetscInt,PetscReal),(F,icntl,val));CHKERRQ(ierr);
14398928b65cSHong Zhang   PetscFunctionReturn(0);
14408928b65cSHong Zhang }
14418928b65cSHong Zhang 
1442bc6112feSHong Zhang #undef __FUNCT__
1443bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetCntl"
1444*a21f80fcSHong Zhang /*@
1445*a21f80fcSHong Zhang   MatMumpsGetCntl - Get MUMPS parameter CNTL()
1446*a21f80fcSHong Zhang 
1447*a21f80fcSHong Zhang    Logically Collective on Mat
1448*a21f80fcSHong Zhang 
1449*a21f80fcSHong Zhang    Input Parameters:
1450*a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1451*a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array CNTL()
1452*a21f80fcSHong Zhang 
1453*a21f80fcSHong Zhang   Output Parameter:
1454*a21f80fcSHong Zhang .  val - value of MUMPS CNTL(icntl)
1455*a21f80fcSHong Zhang 
1456*a21f80fcSHong Zhang    Level: beginner
1457*a21f80fcSHong Zhang 
1458*a21f80fcSHong Zhang    References: MUMPS Users' Guide
1459*a21f80fcSHong Zhang 
1460*a21f80fcSHong Zhang .seealso: MatGetFactor()
1461*a21f80fcSHong Zhang @*/
1462bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl(Mat F,PetscInt icntl,PetscReal *val)
1463bc6112feSHong Zhang {
1464bc6112feSHong Zhang   PetscErrorCode ierr;
1465bc6112feSHong Zhang 
1466bc6112feSHong Zhang   PetscFunctionBegin;
1467bc6112feSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
1468bc6112feSHong Zhang   PetscValidRealPointer(val,3);
1469bc6112feSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetCntl_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
1470bc6112feSHong Zhang   PetscFunctionReturn(0);
1471bc6112feSHong Zhang }
1472bc6112feSHong Zhang 
1473bc6112feSHong Zhang #undef __FUNCT__
1474ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfo_MUMPS"
1475ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F,PetscInt icntl,PetscInt *info)
1476bc6112feSHong Zhang {
1477bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1478bc6112feSHong Zhang 
1479bc6112feSHong Zhang   PetscFunctionBegin;
1480bc6112feSHong Zhang   *info = mumps->id.INFO(icntl);
1481bc6112feSHong Zhang   PetscFunctionReturn(0);
1482bc6112feSHong Zhang }
1483bc6112feSHong Zhang 
1484bc6112feSHong Zhang #undef __FUNCT__
1485ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfog_MUMPS"
1486ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F,PetscInt icntl,PetscInt *infog)
1487bc6112feSHong Zhang {
1488bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1489bc6112feSHong Zhang 
1490bc6112feSHong Zhang   PetscFunctionBegin;
1491bc6112feSHong Zhang   *infog = mumps->id.INFOG(icntl);
1492bc6112feSHong Zhang   PetscFunctionReturn(0);
1493bc6112feSHong Zhang }
1494bc6112feSHong Zhang 
1495bc6112feSHong Zhang #undef __FUNCT__
1496ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfo_MUMPS"
1497ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfo)
1498bc6112feSHong Zhang {
1499bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1500bc6112feSHong Zhang 
1501bc6112feSHong Zhang   PetscFunctionBegin;
1502bc6112feSHong Zhang   *rinfo = mumps->id.RINFO(icntl);
1503bc6112feSHong Zhang   PetscFunctionReturn(0);
1504bc6112feSHong Zhang }
1505bc6112feSHong Zhang 
1506bc6112feSHong Zhang #undef __FUNCT__
1507ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfog_MUMPS"
1508ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfog)
1509bc6112feSHong Zhang {
1510bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1511bc6112feSHong Zhang 
1512bc6112feSHong Zhang   PetscFunctionBegin;
1513bc6112feSHong Zhang   *rinfog = mumps->id.RINFOG(icntl);
1514bc6112feSHong Zhang   PetscFunctionReturn(0);
1515bc6112feSHong Zhang }
1516bc6112feSHong Zhang 
1517bc6112feSHong Zhang #undef __FUNCT__
1518ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfo"
1519*a21f80fcSHong Zhang /*@
1520*a21f80fcSHong Zhang   MatMumpsGetInfo - Get MUMPS parameter INFO()
1521*a21f80fcSHong Zhang 
1522*a21f80fcSHong Zhang    Logically Collective on Mat
1523*a21f80fcSHong Zhang 
1524*a21f80fcSHong Zhang    Input Parameters:
1525*a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1526*a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array INFO()
1527*a21f80fcSHong Zhang 
1528*a21f80fcSHong Zhang   Output Parameter:
1529*a21f80fcSHong Zhang .  ival - value of MUMPS INFO(icntl)
1530*a21f80fcSHong Zhang 
1531*a21f80fcSHong Zhang    Level: beginner
1532*a21f80fcSHong Zhang 
1533*a21f80fcSHong Zhang    References: MUMPS Users' Guide
1534*a21f80fcSHong Zhang 
1535*a21f80fcSHong Zhang .seealso: MatGetFactor()
1536*a21f80fcSHong Zhang @*/
1537ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo(Mat F,PetscInt icntl,PetscInt *ival)
1538bc6112feSHong Zhang {
1539bc6112feSHong Zhang   PetscErrorCode ierr;
1540bc6112feSHong Zhang 
1541bc6112feSHong Zhang   PetscFunctionBegin;
1542ca810319SHong Zhang   PetscValidIntPointer(ival,3);
1543ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetInfo_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
1544bc6112feSHong Zhang   PetscFunctionReturn(0);
1545bc6112feSHong Zhang }
1546bc6112feSHong Zhang 
1547bc6112feSHong Zhang #undef __FUNCT__
1548ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfog"
1549*a21f80fcSHong Zhang /*@
1550*a21f80fcSHong Zhang   MatMumpsGetInfog - Get MUMPS parameter INFOG()
1551*a21f80fcSHong Zhang 
1552*a21f80fcSHong Zhang    Logically Collective on Mat
1553*a21f80fcSHong Zhang 
1554*a21f80fcSHong Zhang    Input Parameters:
1555*a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1556*a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array INFOG()
1557*a21f80fcSHong Zhang 
1558*a21f80fcSHong Zhang   Output Parameter:
1559*a21f80fcSHong Zhang .  ival - value of MUMPS INFOG(icntl)
1560*a21f80fcSHong Zhang 
1561*a21f80fcSHong Zhang    Level: beginner
1562*a21f80fcSHong Zhang 
1563*a21f80fcSHong Zhang    References: MUMPS Users' Guide
1564*a21f80fcSHong Zhang 
1565*a21f80fcSHong Zhang .seealso: MatGetFactor()
1566*a21f80fcSHong Zhang @*/
1567ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog(Mat F,PetscInt icntl,PetscInt *ival)
1568bc6112feSHong Zhang {
1569bc6112feSHong Zhang   PetscErrorCode ierr;
1570bc6112feSHong Zhang 
1571bc6112feSHong Zhang   PetscFunctionBegin;
1572ca810319SHong Zhang   PetscValidIntPointer(ival,3);
1573ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetInfog_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
1574bc6112feSHong Zhang   PetscFunctionReturn(0);
1575bc6112feSHong Zhang }
1576bc6112feSHong Zhang 
1577bc6112feSHong Zhang #undef __FUNCT__
1578ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfo"
1579*a21f80fcSHong Zhang /*@
1580*a21f80fcSHong Zhang   MatMumpsGetRinfo - Get MUMPS parameter RINFO()
1581*a21f80fcSHong Zhang 
1582*a21f80fcSHong Zhang    Logically Collective on Mat
1583*a21f80fcSHong Zhang 
1584*a21f80fcSHong Zhang    Input Parameters:
1585*a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1586*a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array RINFO()
1587*a21f80fcSHong Zhang 
1588*a21f80fcSHong Zhang   Output Parameter:
1589*a21f80fcSHong Zhang .  val - value of MUMPS RINFO(icntl)
1590*a21f80fcSHong Zhang 
1591*a21f80fcSHong Zhang    Level: beginner
1592*a21f80fcSHong Zhang 
1593*a21f80fcSHong Zhang    References: MUMPS Users' Guide
1594*a21f80fcSHong Zhang 
1595*a21f80fcSHong Zhang .seealso: MatGetFactor()
1596*a21f80fcSHong Zhang @*/
1597ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo(Mat F,PetscInt icntl,PetscReal *val)
1598bc6112feSHong Zhang {
1599bc6112feSHong Zhang   PetscErrorCode ierr;
1600bc6112feSHong Zhang 
1601bc6112feSHong Zhang   PetscFunctionBegin;
1602bc6112feSHong Zhang   PetscValidRealPointer(val,3);
1603ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetRinfo_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
1604bc6112feSHong Zhang   PetscFunctionReturn(0);
1605bc6112feSHong Zhang }
1606bc6112feSHong Zhang 
1607bc6112feSHong Zhang #undef __FUNCT__
1608ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfog"
1609*a21f80fcSHong Zhang /*@
1610*a21f80fcSHong Zhang   MatMumpsGetRinfog - Get MUMPS parameter RINFOG()
1611*a21f80fcSHong Zhang 
1612*a21f80fcSHong Zhang    Logically Collective on Mat
1613*a21f80fcSHong Zhang 
1614*a21f80fcSHong Zhang    Input Parameters:
1615*a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1616*a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array RINFOG()
1617*a21f80fcSHong Zhang 
1618*a21f80fcSHong Zhang   Output Parameter:
1619*a21f80fcSHong Zhang .  val - value of MUMPS RINFOG(icntl)
1620*a21f80fcSHong Zhang 
1621*a21f80fcSHong Zhang    Level: beginner
1622*a21f80fcSHong Zhang 
1623*a21f80fcSHong Zhang    References: MUMPS Users' Guide
1624*a21f80fcSHong Zhang 
1625*a21f80fcSHong Zhang .seealso: MatGetFactor()
1626*a21f80fcSHong Zhang @*/
1627ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog(Mat F,PetscInt icntl,PetscReal *val)
1628bc6112feSHong Zhang {
1629bc6112feSHong Zhang   PetscErrorCode ierr;
1630bc6112feSHong Zhang 
1631bc6112feSHong Zhang   PetscFunctionBegin;
1632bc6112feSHong Zhang   PetscValidRealPointer(val,3);
1633ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetRinfog_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
1634bc6112feSHong Zhang   PetscFunctionReturn(0);
1635bc6112feSHong Zhang }
1636bc6112feSHong Zhang 
163724b6179bSKris Buschelman /*MC
16382692d6eeSBarry Smith   MATSOLVERMUMPS -  A matrix type providing direct solvers (LU and Cholesky) for
163924b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
164024b6179bSKris Buschelman 
164141c8de11SBarry Smith   Works with MATAIJ and MATSBAIJ matrices
164224b6179bSKris Buschelman 
164324b6179bSKris Buschelman   Options Database Keys:
1644fb8376fbSHong Zhang + -mat_mumps_icntl_4 <0,...,4> - print level
164524b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
164664e6c443SBarry Smith . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guidec)
164724b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
164824b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
164994b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
165024b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
165124b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
165224b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
165324b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
165424b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
165524b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
165624b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
165724b6179bSKris Buschelman 
165824b6179bSKris Buschelman   Level: beginner
165924b6179bSKris Buschelman 
166041c8de11SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage
166141c8de11SBarry Smith 
166224b6179bSKris Buschelman M*/
166324b6179bSKris Buschelman 
166435bd34faSBarry Smith #undef __FUNCT__
166535bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_mumps"
1666f7a08781SBarry Smith static PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type)
166735bd34faSBarry Smith {
166835bd34faSBarry Smith   PetscFunctionBegin;
16692692d6eeSBarry Smith   *type = MATSOLVERMUMPS;
167035bd34faSBarry Smith   PetscFunctionReturn(0);
167135bd34faSBarry Smith }
167235bd34faSBarry Smith 
1673bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI AIJ matrices */
16742877fffaSHong Zhang #undef __FUNCT__
1675bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_aij_mumps"
16768cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
16772877fffaSHong Zhang {
16782877fffaSHong Zhang   Mat            B;
16792877fffaSHong Zhang   PetscErrorCode ierr;
16802877fffaSHong Zhang   Mat_MUMPS      *mumps;
1681ace3abfcSBarry Smith   PetscBool      isSeqAIJ;
16822877fffaSHong Zhang 
16832877fffaSHong Zhang   PetscFunctionBegin;
16842877fffaSHong Zhang   /* Create the factorization matrix */
1685251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr);
1686ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
16872877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
16882877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1689bccb9932SShri Abhyankar   if (isSeqAIJ) {
16900298fd71SBarry Smith     ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr);
1691bccb9932SShri Abhyankar   } else {
16920298fd71SBarry Smith     ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr);
1693bccb9932SShri Abhyankar   }
16942877fffaSHong Zhang 
1695b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
16962205254eSKarl Rupp 
16972877fffaSHong Zhang   B->ops->view    = MatView_MUMPS;
169835bd34faSBarry Smith   B->ops->getinfo = MatGetInfo_MUMPS;
16992205254eSKarl Rupp 
1700bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
1701bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
1702bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
1703bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
1704bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
1705bc6112feSHong Zhang 
1706ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
1707ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
1708ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
1709ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
1710450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
1711450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
1712d5f3da31SBarry Smith     B->factortype            = MAT_FACTOR_LU;
1713bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
1714bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
1715746480a1SHong Zhang     mumps->sym = 0;
1716dcd589f8SShri Abhyankar   } else {
171767877ebaSShri Abhyankar     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
1718450b117fSShri Abhyankar     B->factortype                  = MAT_FACTOR_CHOLESKY;
1719bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
1720bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
17216fdc2a6dSBarry Smith     if (A->spd_set && A->spd) mumps->sym = 1;
17226fdc2a6dSBarry Smith     else                      mumps->sym = 2;
1723450b117fSShri Abhyankar   }
17242877fffaSHong Zhang 
17252877fffaSHong Zhang   mumps->isAIJ    = PETSC_TRUE;
1726bf0cc555SLisandro Dalcin   mumps->Destroy  = B->ops->destroy;
17272877fffaSHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
17282877fffaSHong Zhang   B->spptr        = (void*)mumps;
17292205254eSKarl Rupp 
1730f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1731746480a1SHong Zhang 
17322877fffaSHong Zhang   *F = B;
17332877fffaSHong Zhang   PetscFunctionReturn(0);
17342877fffaSHong Zhang }
17352877fffaSHong Zhang 
1736bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI SBAIJ matrices */
17372877fffaSHong Zhang #undef __FUNCT__
1738bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_sbaij_mumps"
17398cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
17402877fffaSHong Zhang {
17412877fffaSHong Zhang   Mat            B;
17422877fffaSHong Zhang   PetscErrorCode ierr;
17432877fffaSHong Zhang   Mat_MUMPS      *mumps;
1744ace3abfcSBarry Smith   PetscBool      isSeqSBAIJ;
17452877fffaSHong Zhang 
17462877fffaSHong Zhang   PetscFunctionBegin;
1747ce94432eSBarry Smith   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
1748ce94432eSBarry 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");
1749251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
17502877fffaSHong Zhang   /* Create the factorization matrix */
1751ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
17522877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
17532877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1754b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
1755bccb9932SShri Abhyankar   if (isSeqSBAIJ) {
17560298fd71SBarry Smith     ierr = MatSeqSBAIJSetPreallocation(B,1,0,NULL);CHKERRQ(ierr);
17572205254eSKarl Rupp 
175816ebf90aSShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
1759dcd589f8SShri Abhyankar   } else {
17600298fd71SBarry Smith     ierr = MatMPISBAIJSetPreallocation(B,1,0,NULL,0,NULL);CHKERRQ(ierr);
17612205254eSKarl Rupp 
1762bccb9932SShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
1763bccb9932SShri Abhyankar   }
1764bccb9932SShri Abhyankar 
176567877ebaSShri Abhyankar   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
1766bccb9932SShri Abhyankar   B->ops->view                   = MatView_MUMPS;
17672205254eSKarl Rupp 
1768bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
1769b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
1770b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
1771b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
1772b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
1773bc6112feSHong Zhang 
1774ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
1775ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
1776ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
1777ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
17782205254eSKarl Rupp 
1779f4762488SHong Zhang   B->factortype = MAT_FACTOR_CHOLESKY;
17806fdc2a6dSBarry Smith   if (A->spd_set && A->spd) mumps->sym = 1;
17816fdc2a6dSBarry Smith   else                      mumps->sym = 2;
1782a214ac2aSShri Abhyankar 
1783bccb9932SShri Abhyankar   mumps->isAIJ    = PETSC_FALSE;
1784bf0cc555SLisandro Dalcin   mumps->Destroy  = B->ops->destroy;
1785f3c0ef26SHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
17862877fffaSHong Zhang   B->spptr        = (void*)mumps;
17872205254eSKarl Rupp 
1788f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1789746480a1SHong Zhang 
17902877fffaSHong Zhang   *F = B;
17912877fffaSHong Zhang   PetscFunctionReturn(0);
17922877fffaSHong Zhang }
179397969023SHong Zhang 
1794450b117fSShri Abhyankar #undef __FUNCT__
1795bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_baij_mumps"
17968cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
179767877ebaSShri Abhyankar {
179867877ebaSShri Abhyankar   Mat            B;
179967877ebaSShri Abhyankar   PetscErrorCode ierr;
180067877ebaSShri Abhyankar   Mat_MUMPS      *mumps;
1801ace3abfcSBarry Smith   PetscBool      isSeqBAIJ;
180267877ebaSShri Abhyankar 
180367877ebaSShri Abhyankar   PetscFunctionBegin;
180467877ebaSShri Abhyankar   /* Create the factorization matrix */
1805251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr);
1806ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
180767877ebaSShri Abhyankar   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
180867877ebaSShri Abhyankar   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1809bccb9932SShri Abhyankar   if (isSeqBAIJ) {
18100298fd71SBarry Smith     ierr = MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,NULL);CHKERRQ(ierr);
1811bccb9932SShri Abhyankar   } else {
18120298fd71SBarry Smith     ierr = MatMPIBAIJSetPreallocation(B,A->rmap->bs,0,NULL,0,NULL);CHKERRQ(ierr);
1813bccb9932SShri Abhyankar   }
1814450b117fSShri Abhyankar 
1815b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
1816450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
1817450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
1818450b117fSShri Abhyankar     B->factortype            = MAT_FACTOR_LU;
1819bccb9932SShri Abhyankar     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
1820bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
1821746480a1SHong Zhang     mumps->sym = 0;
1822f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n");
1823bccb9932SShri Abhyankar 
1824450b117fSShri Abhyankar   B->ops->view = MatView_MUMPS;
18252205254eSKarl Rupp 
1826bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
1827bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
1828bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
1829bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
1830bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
1831bc6112feSHong Zhang 
1832ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
1833ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
1834ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
1835ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
1836450b117fSShri Abhyankar 
1837450b117fSShri Abhyankar   mumps->isAIJ    = PETSC_TRUE;
1838bf0cc555SLisandro Dalcin   mumps->Destroy  = B->ops->destroy;
1839450b117fSShri Abhyankar   B->ops->destroy = MatDestroy_MUMPS;
1840450b117fSShri Abhyankar   B->spptr        = (void*)mumps;
18412205254eSKarl Rupp 
1842f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1843746480a1SHong Zhang 
1844450b117fSShri Abhyankar   *F = B;
1845450b117fSShri Abhyankar   PetscFunctionReturn(0);
1846450b117fSShri Abhyankar }
1847