xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision 2907cef92103ca874bcf8ea02a2fe5afb9ad785d)
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)
11*2907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
12*2907cef9SHong Zhang #include <cmumps_c.h>
13*2907cef9SHong Zhang #else
14c6db04a5SJed Brown #include <zmumps_c.h>
15*2907cef9SHong Zhang #endif
16*2907cef9SHong Zhang #else
17*2907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
18*2907cef9SHong Zhang #include <smumps_c.h>
19397b6df1SKris Buschelman #else
20c6db04a5SJed Brown #include <dmumps_c.h>
21397b6df1SKris Buschelman #endif
22*2907cef9SHong 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 
30*2907cef9SHong Zhang /* calls to MUMPS */
31*2907cef9SHong Zhang #if defined(PETSC_USE_COMPLEX)
32*2907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
33*2907cef9SHong Zhang #define PetscMUMPS_c cmumps_c
34*2907cef9SHong Zhang #else
35*2907cef9SHong Zhang #define PetscMUMPS_c zmumps_c
36*2907cef9SHong Zhang #endif
37*2907cef9SHong Zhang #else
38*2907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
39*2907cef9SHong Zhang #define PetscMUMPS_c smumps_c
40*2907cef9SHong Zhang #else
41*2907cef9SHong Zhang #define PetscMUMPS_c dmumps_c
42*2907cef9SHong Zhang #endif
43*2907cef9SHong Zhang #endif
44*2907cef9SHong 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)
56*2907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
57*2907cef9SHong Zhang   CMUMPS_STRUC_C id;
58*2907cef9SHong Zhang #else
59397b6df1SKris Buschelman   ZMUMPS_STRUC_C id;
60*2907cef9SHong Zhang #endif
61*2907cef9SHong Zhang #else
62*2907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
63*2907cef9SHong Zhang   SMUMPS_STRUC_C id;
64397b6df1SKris Buschelman #else
65397b6df1SKris Buschelman   DMUMPS_STRUC_C id;
66397b6df1SKris Buschelman #endif
67*2907cef9SHong Zhang #endif
68*2907cef9SHong Zhang 
69397b6df1SKris Buschelman   MatStructure   matstruc;
70c1490034SHong Zhang   PetscMPIInt    myid,size;
7116ebf90aSShri Abhyankar   PetscInt       *irn,*jcn,nz,sym,nSolve;
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;
77bf0cc555SLisandro Dalcin   PetscErrorCode (*Destroy)(Mat);
78bccb9932SShri Abhyankar   PetscErrorCode (*ConvertToTriples)(Mat, int, MatReuse, int*, int**, int**, PetscScalar**);
79f0c56d0fSKris Buschelman } Mat_MUMPS;
80f0c56d0fSKris Buschelman 
8109573ac7SBarry Smith extern PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);
82b24902e0SBarry Smith 
8367877ebaSShri Abhyankar 
8467877ebaSShri Abhyankar /* MatConvertToTriples_A_B */
8567877ebaSShri Abhyankar /*convert Petsc matrix to triples: row[nz], col[nz], val[nz] */
86397b6df1SKris Buschelman /*
87397b6df1SKris Buschelman   input:
8867877ebaSShri Abhyankar     A       - matrix in aij,baij or sbaij (bs=1) format
89397b6df1SKris Buschelman     shift   - 0: C style output triple; 1: Fortran style output triple.
90bccb9932SShri Abhyankar     reuse   - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple
91bccb9932SShri Abhyankar               MAT_REUSE_MATRIX:   only the values in v array are updated
92397b6df1SKris Buschelman   output:
93397b6df1SKris Buschelman     nnz     - dim of r, c, and v (number of local nonzero entries of A)
94397b6df1SKris Buschelman     r, c, v - row and col index, matrix values (matrix triples)
95397b6df1SKris Buschelman  */
9616ebf90aSShri Abhyankar 
9716ebf90aSShri Abhyankar #undef __FUNCT__
9816ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqaij"
99bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
100b24902e0SBarry Smith {
101185f6596SHong Zhang   const PetscInt   *ai,*aj,*ajj,M=A->rmap->n;
10267877ebaSShri Abhyankar   PetscInt         nz,rnz,i,j;
103dfbe8321SBarry Smith   PetscErrorCode   ierr;
104c1490034SHong Zhang   PetscInt         *row,*col;
10516ebf90aSShri Abhyankar   Mat_SeqAIJ       *aa=(Mat_SeqAIJ*)A->data;
106397b6df1SKris Buschelman 
107397b6df1SKris Buschelman   PetscFunctionBegin;
10816ebf90aSShri Abhyankar   *v=aa->a;
109bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX){
11016ebf90aSShri Abhyankar     nz = aa->nz; ai = aa->i; aj = aa->j;
11116ebf90aSShri Abhyankar     *nnz = nz;
112185f6596SHong Zhang     ierr = PetscMalloc(2*nz*sizeof(PetscInt), &row);CHKERRQ(ierr);
113185f6596SHong Zhang     col  = row + nz;
114185f6596SHong Zhang 
11516ebf90aSShri Abhyankar     nz = 0;
11616ebf90aSShri Abhyankar     for (i=0; i<M; i++) {
11716ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
11867877ebaSShri Abhyankar       ajj = aj + ai[i];
11967877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
12067877ebaSShri Abhyankar 	row[nz] = i+shift; col[nz++] = ajj[j] + shift;
12116ebf90aSShri Abhyankar       }
12216ebf90aSShri Abhyankar     }
12316ebf90aSShri Abhyankar     *r = row; *c = col;
12416ebf90aSShri Abhyankar   }
12516ebf90aSShri Abhyankar   PetscFunctionReturn(0);
12616ebf90aSShri Abhyankar }
127397b6df1SKris Buschelman 
12816ebf90aSShri Abhyankar #undef __FUNCT__
12967877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqbaij_seqaij"
130bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
13167877ebaSShri Abhyankar {
13267877ebaSShri Abhyankar   Mat_SeqBAIJ        *aa=(Mat_SeqBAIJ*)A->data;
13367877ebaSShri Abhyankar   const PetscInt     *ai,*aj,*ajj,bs=A->rmap->bs,bs2=aa->bs2,M=A->rmap->N/bs;
1340ad0caddSJed Brown   PetscInt           nz,idx=0,rnz,i,j,k,m;
13567877ebaSShri Abhyankar   PetscErrorCode     ierr;
13667877ebaSShri Abhyankar   PetscInt           *row,*col;
13767877ebaSShri Abhyankar 
13867877ebaSShri Abhyankar   PetscFunctionBegin;
139cf3759fdSShri Abhyankar   *v = aa->a;
140bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX){
141cf3759fdSShri Abhyankar     ai = aa->i; aj = aa->j;
14267877ebaSShri Abhyankar     nz = bs2*aa->nz;
14367877ebaSShri Abhyankar     *nnz = nz;
144185f6596SHong Zhang     ierr = PetscMalloc(2*nz*sizeof(PetscInt), &row);CHKERRQ(ierr);
145185f6596SHong Zhang     col  = row + nz;
146185f6596SHong Zhang 
14767877ebaSShri Abhyankar     for (i=0; i<M; i++) {
14867877ebaSShri Abhyankar       ajj = aj + ai[i];
14967877ebaSShri Abhyankar       rnz = ai[i+1] - ai[i];
15067877ebaSShri Abhyankar       for (k=0; k<rnz; k++) {
15167877ebaSShri Abhyankar 	for (j=0; j<bs; j++) {
15267877ebaSShri Abhyankar 	  for (m=0; m<bs; m++) {
15367877ebaSShri Abhyankar 	    row[idx]     = i*bs + m + shift;
154cf3759fdSShri Abhyankar 	    col[idx++]   = bs*(ajj[k]) + j + shift;
15567877ebaSShri Abhyankar 	  }
15667877ebaSShri Abhyankar 	}
15767877ebaSShri Abhyankar       }
15867877ebaSShri Abhyankar     }
159cf3759fdSShri Abhyankar     *r = row; *c = col;
16067877ebaSShri Abhyankar   }
16167877ebaSShri Abhyankar   PetscFunctionReturn(0);
16267877ebaSShri Abhyankar }
16367877ebaSShri Abhyankar 
16467877ebaSShri Abhyankar #undef __FUNCT__
16516ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqsbaij_seqsbaij"
166bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
16716ebf90aSShri Abhyankar {
16867877ebaSShri Abhyankar   const PetscInt   *ai, *aj,*ajj,M=A->rmap->n;
16967877ebaSShri Abhyankar   PetscInt         nz,rnz,i,j;
17016ebf90aSShri Abhyankar   PetscErrorCode   ierr;
17116ebf90aSShri Abhyankar   PetscInt         *row,*col;
17216ebf90aSShri Abhyankar   Mat_SeqSBAIJ     *aa=(Mat_SeqSBAIJ*)A->data;
17316ebf90aSShri Abhyankar 
17416ebf90aSShri Abhyankar   PetscFunctionBegin;
175882afa5aSHong Zhang   *v = aa->a;
176bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX){
17716ebf90aSShri Abhyankar     nz = aa->nz;ai=aa->i; aj=aa->j;*v=aa->a;
17816ebf90aSShri Abhyankar     *nnz = nz;
179185f6596SHong Zhang     ierr = PetscMalloc(2*nz*sizeof(PetscInt), &row);CHKERRQ(ierr);
180185f6596SHong Zhang     col  = row + nz;
181185f6596SHong Zhang 
18216ebf90aSShri Abhyankar     nz = 0;
18316ebf90aSShri Abhyankar     for (i=0; i<M; i++) {
18416ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
18567877ebaSShri Abhyankar       ajj = aj + ai[i];
18667877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
18767877ebaSShri Abhyankar 	row[nz] = i+shift; col[nz++] = ajj[j] + shift;
18816ebf90aSShri Abhyankar       }
18916ebf90aSShri Abhyankar     }
19016ebf90aSShri Abhyankar     *r = row; *c = col;
19116ebf90aSShri Abhyankar   }
19216ebf90aSShri Abhyankar   PetscFunctionReturn(0);
19316ebf90aSShri Abhyankar }
19416ebf90aSShri Abhyankar 
19516ebf90aSShri Abhyankar #undef __FUNCT__
19616ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqsbaij"
197bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
19816ebf90aSShri Abhyankar {
19967877ebaSShri Abhyankar   const PetscInt     *ai,*aj,*ajj,*adiag,M=A->rmap->n;
20067877ebaSShri Abhyankar   PetscInt           nz,rnz,i,j;
20167877ebaSShri Abhyankar   const PetscScalar  *av,*v1;
20216ebf90aSShri Abhyankar   PetscScalar        *val;
20316ebf90aSShri Abhyankar   PetscErrorCode     ierr;
20416ebf90aSShri Abhyankar   PetscInt           *row,*col;
20516ebf90aSShri Abhyankar   Mat_SeqSBAIJ       *aa=(Mat_SeqSBAIJ*)A->data;
20616ebf90aSShri Abhyankar 
20716ebf90aSShri Abhyankar   PetscFunctionBegin;
20816ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j;av=aa->a;
20916ebf90aSShri Abhyankar   adiag=aa->diag;
210bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX){
21116ebf90aSShri Abhyankar     nz = M + (aa->nz-M)/2;
21216ebf90aSShri Abhyankar     *nnz = nz;
213185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
214185f6596SHong Zhang     col  = row + nz;
215185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
216185f6596SHong Zhang 
21716ebf90aSShri Abhyankar     nz = 0;
21816ebf90aSShri Abhyankar     for (i=0; i<M; i++) {
21916ebf90aSShri Abhyankar       rnz = ai[i+1] - adiag[i];
22067877ebaSShri Abhyankar       ajj  = aj + adiag[i];
221cf3759fdSShri Abhyankar       v1   = av + adiag[i];
22267877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
22367877ebaSShri Abhyankar 	row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
22416ebf90aSShri Abhyankar       }
22516ebf90aSShri Abhyankar     }
22616ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
227397b6df1SKris Buschelman   } else {
22816ebf90aSShri Abhyankar     nz = 0; val = *v;
22916ebf90aSShri Abhyankar     for (i=0; i <M; i++) {
23016ebf90aSShri Abhyankar       rnz = ai[i+1] - adiag[i];
23167877ebaSShri Abhyankar       ajj = aj + adiag[i];
23267877ebaSShri Abhyankar       v1  = av + adiag[i];
23367877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
23467877ebaSShri Abhyankar 	val[nz++] = v1[j];
23516ebf90aSShri Abhyankar       }
23616ebf90aSShri Abhyankar     }
23716ebf90aSShri Abhyankar   }
23816ebf90aSShri Abhyankar   PetscFunctionReturn(0);
23916ebf90aSShri Abhyankar }
24016ebf90aSShri Abhyankar 
24116ebf90aSShri Abhyankar #undef __FUNCT__
24216ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpisbaij_mpisbaij"
243bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
24416ebf90aSShri Abhyankar {
24516ebf90aSShri Abhyankar   const PetscInt     *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
24616ebf90aSShri Abhyankar   PetscErrorCode     ierr;
24716ebf90aSShri Abhyankar   PetscInt           rstart,nz,i,j,jj,irow,countA,countB;
24816ebf90aSShri Abhyankar   PetscInt           *row,*col;
24916ebf90aSShri Abhyankar   const PetscScalar  *av, *bv,*v1,*v2;
25016ebf90aSShri Abhyankar   PetscScalar        *val;
251397b6df1SKris Buschelman   Mat_MPISBAIJ       *mat =  (Mat_MPISBAIJ*)A->data;
252397b6df1SKris Buschelman   Mat_SeqSBAIJ       *aa=(Mat_SeqSBAIJ*)(mat->A)->data;
253397b6df1SKris Buschelman   Mat_SeqBAIJ        *bb=(Mat_SeqBAIJ*)(mat->B)->data;
25416ebf90aSShri Abhyankar 
25516ebf90aSShri Abhyankar   PetscFunctionBegin;
256d0f46423SBarry Smith   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
257397b6df1SKris Buschelman   garray = mat->garray;
258397b6df1SKris Buschelman   av=aa->a; bv=bb->a;
259397b6df1SKris Buschelman 
260bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX){
26116ebf90aSShri Abhyankar     nz = aa->nz + bb->nz;
26216ebf90aSShri Abhyankar     *nnz = nz;
263185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
264185f6596SHong Zhang     col  = row + nz;
265185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
266185f6596SHong Zhang 
267397b6df1SKris Buschelman     *r = row; *c = col; *v = val;
268397b6df1SKris Buschelman   } else {
269397b6df1SKris Buschelman     row = *r; col = *c; val = *v;
270397b6df1SKris Buschelman   }
271397b6df1SKris Buschelman 
272028e57e8SHong Zhang   jj = 0; irow = rstart;
273397b6df1SKris Buschelman   for ( i=0; i<m; i++ ) {
274397b6df1SKris Buschelman     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
275397b6df1SKris Buschelman     countA = ai[i+1] - ai[i];
276397b6df1SKris Buschelman     countB = bi[i+1] - bi[i];
277397b6df1SKris Buschelman     bjj    = bj + bi[i];
27816ebf90aSShri Abhyankar     v1     = av + ai[i];
27916ebf90aSShri Abhyankar     v2     = bv + bi[i];
280397b6df1SKris Buschelman 
281397b6df1SKris Buschelman     /* A-part */
282397b6df1SKris Buschelman     for (j=0; j<countA; j++){
283bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
284397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
285397b6df1SKris Buschelman       }
28616ebf90aSShri Abhyankar       val[jj++] = v1[j];
287397b6df1SKris Buschelman     }
28816ebf90aSShri Abhyankar 
28916ebf90aSShri Abhyankar     /* B-part */
29016ebf90aSShri Abhyankar     for (j=0; j < countB; j++){
291bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
292397b6df1SKris Buschelman 	row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
293397b6df1SKris Buschelman       }
29416ebf90aSShri Abhyankar       val[jj++] = v2[j];
29516ebf90aSShri Abhyankar     }
29616ebf90aSShri Abhyankar     irow++;
29716ebf90aSShri Abhyankar   }
29816ebf90aSShri Abhyankar   PetscFunctionReturn(0);
29916ebf90aSShri Abhyankar }
30016ebf90aSShri Abhyankar 
30116ebf90aSShri Abhyankar #undef __FUNCT__
30216ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpiaij"
303bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
30416ebf90aSShri Abhyankar {
30516ebf90aSShri Abhyankar   const PetscInt     *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
30616ebf90aSShri Abhyankar   PetscErrorCode     ierr;
30716ebf90aSShri Abhyankar   PetscInt           rstart,nz,i,j,jj,irow,countA,countB;
30816ebf90aSShri Abhyankar   PetscInt           *row,*col;
30916ebf90aSShri Abhyankar   const PetscScalar  *av, *bv,*v1,*v2;
31016ebf90aSShri Abhyankar   PetscScalar        *val;
31116ebf90aSShri Abhyankar   Mat_MPIAIJ         *mat =  (Mat_MPIAIJ*)A->data;
31216ebf90aSShri Abhyankar   Mat_SeqAIJ         *aa=(Mat_SeqAIJ*)(mat->A)->data;
31316ebf90aSShri Abhyankar   Mat_SeqAIJ         *bb=(Mat_SeqAIJ*)(mat->B)->data;
31416ebf90aSShri Abhyankar 
31516ebf90aSShri Abhyankar   PetscFunctionBegin;
31616ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
31716ebf90aSShri Abhyankar   garray = mat->garray;
31816ebf90aSShri Abhyankar   av=aa->a; bv=bb->a;
31916ebf90aSShri Abhyankar 
320bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX){
32116ebf90aSShri Abhyankar     nz = aa->nz + bb->nz;
32216ebf90aSShri Abhyankar     *nnz = nz;
323185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
324185f6596SHong Zhang     col  = row + nz;
325185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
326185f6596SHong Zhang 
32716ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
32816ebf90aSShri Abhyankar   } else {
32916ebf90aSShri Abhyankar     row = *r; col = *c; val = *v;
33016ebf90aSShri Abhyankar   }
33116ebf90aSShri Abhyankar 
33216ebf90aSShri Abhyankar   jj = 0; irow = rstart;
33316ebf90aSShri Abhyankar   for ( i=0; i<m; i++ ) {
33416ebf90aSShri Abhyankar     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
33516ebf90aSShri Abhyankar     countA = ai[i+1] - ai[i];
33616ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
33716ebf90aSShri Abhyankar     bjj    = bj + bi[i];
33816ebf90aSShri Abhyankar     v1     = av + ai[i];
33916ebf90aSShri Abhyankar     v2     = bv + bi[i];
34016ebf90aSShri Abhyankar 
34116ebf90aSShri Abhyankar     /* A-part */
34216ebf90aSShri Abhyankar     for (j=0; j<countA; j++){
343bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX){
34416ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
34516ebf90aSShri Abhyankar       }
34616ebf90aSShri Abhyankar       val[jj++] = v1[j];
34716ebf90aSShri Abhyankar     }
34816ebf90aSShri Abhyankar 
34916ebf90aSShri Abhyankar     /* B-part */
35016ebf90aSShri Abhyankar     for (j=0; j < countB; j++){
351bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX){
35216ebf90aSShri Abhyankar 	row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
35316ebf90aSShri Abhyankar       }
35416ebf90aSShri Abhyankar       val[jj++] = v2[j];
35516ebf90aSShri Abhyankar     }
35616ebf90aSShri Abhyankar     irow++;
35716ebf90aSShri Abhyankar   }
35816ebf90aSShri Abhyankar   PetscFunctionReturn(0);
35916ebf90aSShri Abhyankar }
36016ebf90aSShri Abhyankar 
36116ebf90aSShri Abhyankar #undef __FUNCT__
36267877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpibaij_mpiaij"
363bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
36467877ebaSShri Abhyankar {
36567877ebaSShri Abhyankar   Mat_MPIBAIJ        *mat =  (Mat_MPIBAIJ*)A->data;
36667877ebaSShri Abhyankar   Mat_SeqBAIJ        *aa=(Mat_SeqBAIJ*)(mat->A)->data;
36767877ebaSShri Abhyankar   Mat_SeqBAIJ        *bb=(Mat_SeqBAIJ*)(mat->B)->data;
36867877ebaSShri Abhyankar   const PetscInt     *ai = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj;
369d985c460SShri Abhyankar   const PetscInt     *garray = mat->garray,mbs=mat->mbs,rstart=A->rmap->rstart;
37067877ebaSShri Abhyankar   const PetscInt     bs = A->rmap->bs,bs2=mat->bs2;
37167877ebaSShri Abhyankar   PetscErrorCode     ierr;
37267877ebaSShri Abhyankar   PetscInt           nz,i,j,k,n,jj,irow,countA,countB,idx;
37367877ebaSShri Abhyankar   PetscInt           *row,*col;
37467877ebaSShri Abhyankar   const PetscScalar  *av=aa->a, *bv=bb->a,*v1,*v2;
37567877ebaSShri Abhyankar   PetscScalar        *val;
37667877ebaSShri Abhyankar 
37767877ebaSShri Abhyankar   PetscFunctionBegin;
37867877ebaSShri Abhyankar 
379bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
38067877ebaSShri Abhyankar     nz = bs2*(aa->nz + bb->nz);
38167877ebaSShri Abhyankar     *nnz = nz;
382185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
383185f6596SHong Zhang     col  = row + nz;
384185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
385185f6596SHong Zhang 
38667877ebaSShri Abhyankar     *r = row; *c = col; *v = val;
38767877ebaSShri Abhyankar   } else {
38867877ebaSShri Abhyankar     row = *r; col = *c; val = *v;
38967877ebaSShri Abhyankar   }
39067877ebaSShri Abhyankar 
391d985c460SShri Abhyankar   jj = 0; irow = rstart;
39267877ebaSShri Abhyankar   for ( i=0; i<mbs; i++ ) {
39367877ebaSShri Abhyankar     countA = ai[i+1] - ai[i];
39467877ebaSShri Abhyankar     countB = bi[i+1] - bi[i];
39567877ebaSShri Abhyankar     ajj    = aj + ai[i];
39667877ebaSShri Abhyankar     bjj    = bj + bi[i];
39767877ebaSShri Abhyankar     v1     = av + bs2*ai[i];
39867877ebaSShri Abhyankar     v2     = bv + bs2*bi[i];
39967877ebaSShri Abhyankar 
40067877ebaSShri Abhyankar     idx = 0;
40167877ebaSShri Abhyankar     /* A-part */
40267877ebaSShri Abhyankar     for (k=0; k<countA; k++){
40367877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
40467877ebaSShri Abhyankar 	for (n=0; n<bs; n++) {
405bccb9932SShri Abhyankar 	  if (reuse == MAT_INITIAL_MATRIX){
406d985c460SShri Abhyankar 	    row[jj] = irow + n + shift;
407d985c460SShri Abhyankar 	    col[jj] = rstart + bs*ajj[k] + j + shift;
40867877ebaSShri Abhyankar 	  }
40967877ebaSShri Abhyankar 	  val[jj++] = v1[idx++];
41067877ebaSShri Abhyankar 	}
41167877ebaSShri Abhyankar       }
41267877ebaSShri Abhyankar     }
41367877ebaSShri Abhyankar 
41467877ebaSShri Abhyankar     idx = 0;
41567877ebaSShri Abhyankar     /* B-part */
41667877ebaSShri Abhyankar     for (k=0; k<countB; k++){
41767877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
41867877ebaSShri Abhyankar 	for (n=0; n<bs; n++) {
419bccb9932SShri Abhyankar 	  if (reuse == MAT_INITIAL_MATRIX){
420d985c460SShri Abhyankar 	    row[jj] = irow + n + shift;
421d985c460SShri Abhyankar 	    col[jj] = bs*garray[bjj[k]] + j + shift;
42267877ebaSShri Abhyankar 	  }
423d985c460SShri Abhyankar 	  val[jj++] = v2[idx++];
42467877ebaSShri Abhyankar 	}
42567877ebaSShri Abhyankar       }
42667877ebaSShri Abhyankar     }
427d985c460SShri Abhyankar     irow += bs;
42867877ebaSShri Abhyankar   }
42967877ebaSShri Abhyankar   PetscFunctionReturn(0);
43067877ebaSShri Abhyankar }
43167877ebaSShri Abhyankar 
43267877ebaSShri Abhyankar #undef __FUNCT__
43316ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpisbaij"
434bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
43516ebf90aSShri Abhyankar {
43616ebf90aSShri Abhyankar   const PetscInt     *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
43716ebf90aSShri Abhyankar   PetscErrorCode     ierr;
438e0bace9bSHong Zhang   PetscInt           rstart,nz,nza,nzb,i,j,jj,irow,countA,countB;
43916ebf90aSShri Abhyankar   PetscInt           *row,*col;
44016ebf90aSShri Abhyankar   const PetscScalar  *av, *bv,*v1,*v2;
44116ebf90aSShri Abhyankar   PetscScalar        *val;
44216ebf90aSShri Abhyankar   Mat_MPIAIJ         *mat =  (Mat_MPIAIJ*)A->data;
44316ebf90aSShri Abhyankar   Mat_SeqAIJ         *aa=(Mat_SeqAIJ*)(mat->A)->data;
44416ebf90aSShri Abhyankar   Mat_SeqAIJ         *bb=(Mat_SeqAIJ*)(mat->B)->data;
44516ebf90aSShri Abhyankar 
44616ebf90aSShri Abhyankar   PetscFunctionBegin;
44716ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j; adiag=aa->diag;
44816ebf90aSShri Abhyankar   bi=bb->i; bj=bb->j; garray = mat->garray;
44916ebf90aSShri Abhyankar   av=aa->a; bv=bb->a;
45016ebf90aSShri Abhyankar   rstart = A->rmap->rstart;
45116ebf90aSShri Abhyankar 
452bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
453e0bace9bSHong Zhang     nza = 0;    /* num of upper triangular entries in mat->A, including diagonals */
454e0bace9bSHong Zhang     nzb = 0;    /* num of upper triangular entries in mat->B */
45516ebf90aSShri Abhyankar     for (i=0; i<m; i++){
456e0bace9bSHong Zhang       nza    += (ai[i+1] - adiag[i]);
45716ebf90aSShri Abhyankar       countB  = bi[i+1] - bi[i];
45816ebf90aSShri Abhyankar       bjj     = bj + bi[i];
459e0bace9bSHong Zhang       for (j=0; j<countB; j++){
460e0bace9bSHong Zhang         if (garray[bjj[j]] > rstart) nzb++;
461e0bace9bSHong Zhang       }
462e0bace9bSHong Zhang     }
46316ebf90aSShri Abhyankar 
464e0bace9bSHong Zhang     nz = nza + nzb; /* total nz of upper triangular part of mat */
46516ebf90aSShri Abhyankar     *nnz = nz;
466185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
467185f6596SHong Zhang     col  = row + nz;
468185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
469185f6596SHong Zhang 
47016ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
47116ebf90aSShri Abhyankar   } else {
47216ebf90aSShri Abhyankar     row = *r; col = *c; val = *v;
47316ebf90aSShri Abhyankar   }
47416ebf90aSShri Abhyankar 
47516ebf90aSShri Abhyankar   jj = 0; irow = rstart;
47616ebf90aSShri Abhyankar   for ( i=0; i<m; i++ ) {
47716ebf90aSShri Abhyankar     ajj    = aj + adiag[i];                 /* ptr to the beginning of the diagonal of this row */
47816ebf90aSShri Abhyankar     v1     = av + adiag[i];
47916ebf90aSShri Abhyankar     countA = ai[i+1] - adiag[i];
48016ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
48116ebf90aSShri Abhyankar     bjj    = bj + bi[i];
48216ebf90aSShri Abhyankar     v2     = bv + bi[i];
48316ebf90aSShri Abhyankar 
48416ebf90aSShri Abhyankar      /* A-part */
48516ebf90aSShri Abhyankar     for (j=0; j<countA; j++){
486bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
48716ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
48816ebf90aSShri Abhyankar       }
48916ebf90aSShri Abhyankar       val[jj++] = v1[j];
49016ebf90aSShri Abhyankar     }
49116ebf90aSShri Abhyankar 
49216ebf90aSShri Abhyankar     /* B-part */
49316ebf90aSShri Abhyankar     for (j=0; j < countB; j++){
49416ebf90aSShri Abhyankar       if (garray[bjj[j]] > rstart) {
495bccb9932SShri Abhyankar 	if (reuse == MAT_INITIAL_MATRIX) {
49616ebf90aSShri Abhyankar 	  row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
49716ebf90aSShri Abhyankar 	}
49816ebf90aSShri Abhyankar 	val[jj++] = v2[j];
49916ebf90aSShri Abhyankar       }
500397b6df1SKris Buschelman     }
501397b6df1SKris Buschelman     irow++;
502397b6df1SKris Buschelman   }
503397b6df1SKris Buschelman   PetscFunctionReturn(0);
504397b6df1SKris Buschelman }
505397b6df1SKris Buschelman 
506397b6df1SKris Buschelman #undef __FUNCT__
5073924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS"
508dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A)
509dfbe8321SBarry Smith {
510f0c56d0fSKris Buschelman   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
511dfbe8321SBarry Smith   PetscErrorCode ierr;
512b24902e0SBarry Smith 
513397b6df1SKris Buschelman   PetscFunctionBegin;
514bf0cc555SLisandro Dalcin   if (lu && lu->CleanUpMUMPS) {
515397b6df1SKris Buschelman     /* Terminate instance, deallocate memories */
5166bf464f9SBarry Smith     ierr = PetscFree2(lu->id.sol_loc,lu->id.isol_loc);CHKERRQ(ierr);
5176bf464f9SBarry Smith     ierr = VecScatterDestroy(&lu->scat_rhs);CHKERRQ(ierr);
5186bf464f9SBarry Smith     ierr = VecDestroy(&lu->b_seq);CHKERRQ(ierr);
519bf0cc555SLisandro Dalcin     ierr = VecScatterDestroy(&lu->scat_sol);CHKERRQ(ierr);
520bf0cc555SLisandro Dalcin     ierr = VecDestroy(&lu->x_seq);CHKERRQ(ierr);
5216bf464f9SBarry Smith     ierr=PetscFree(lu->id.perm_in);CHKERRQ(ierr);
522185f6596SHong Zhang     ierr = PetscFree(lu->irn);CHKERRQ(ierr);
523397b6df1SKris Buschelman     lu->id.job=JOB_END;
524*2907cef9SHong Zhang     PetscMUMPS_c(&lu->id);
525397b6df1SKris Buschelman     ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr);
526397b6df1SKris Buschelman   }
527bf0cc555SLisandro Dalcin   if (lu && lu->Destroy) {
528bf0cc555SLisandro Dalcin     ierr = (lu->Destroy)(A);CHKERRQ(ierr);
529bf0cc555SLisandro Dalcin   }
530bf0cc555SLisandro Dalcin   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
531bf0cc555SLisandro Dalcin 
53297969023SHong Zhang   /* clear composed functions */
53397969023SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatFactorGetSolverPackage_C","",PETSC_NULL);CHKERRQ(ierr);
534f250808bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMumpsSetIcntl_C","",PETSC_NULL);CHKERRQ(ierr);
535397b6df1SKris Buschelman   PetscFunctionReturn(0);
536397b6df1SKris Buschelman }
537397b6df1SKris Buschelman 
538397b6df1SKris Buschelman #undef __FUNCT__
539f6c57405SHong Zhang #define __FUNCT__ "MatSolve_MUMPS"
540b24902e0SBarry Smith PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x)
541b24902e0SBarry Smith {
542f0c56d0fSKris Buschelman   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
543d54de34fSKris Buschelman   PetscScalar    *array;
54467877ebaSShri Abhyankar   Vec            b_seq;
545329ec9b3SHong Zhang   IS             is_iden,is_petsc;
546dfbe8321SBarry Smith   PetscErrorCode ierr;
547329ec9b3SHong Zhang   PetscInt       i;
548397b6df1SKris Buschelman 
549397b6df1SKris Buschelman   PetscFunctionBegin;
550329ec9b3SHong Zhang   lu->id.nrhs = 1;
55167877ebaSShri Abhyankar   b_seq = lu->b_seq;
552397b6df1SKris Buschelman   if (lu->size > 1){
553329ec9b3SHong Zhang     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
55467877ebaSShri Abhyankar     ierr = VecScatterBegin(lu->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
55567877ebaSShri Abhyankar     ierr = VecScatterEnd(lu->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
55667877ebaSShri Abhyankar     if (!lu->myid) {ierr = VecGetArray(b_seq,&array);CHKERRQ(ierr);}
557397b6df1SKris Buschelman   } else {  /* size == 1 */
558397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
559397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
560397b6df1SKris Buschelman   }
561397b6df1SKris Buschelman   if (!lu->myid) { /* define rhs on the host */
5628278f211SHong Zhang     lu->id.nrhs = 1;
563397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
564*2907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
565*2907cef9SHong Zhang     lu->id.rhs = (mumps_complex*)array;
566*2907cef9SHong Zhang #else
567397b6df1SKris Buschelman     lu->id.rhs = (mumps_double_complex*)array;
568*2907cef9SHong Zhang #endif
569397b6df1SKris Buschelman #else
570397b6df1SKris Buschelman     lu->id.rhs = array;
571397b6df1SKris Buschelman #endif
572397b6df1SKris Buschelman   }
573397b6df1SKris Buschelman 
574397b6df1SKris Buschelman   /* solve phase */
575329ec9b3SHong Zhang   /*-------------*/
5763d472b54SHong Zhang   lu->id.job = JOB_SOLVE;
577*2907cef9SHong Zhang   PetscMUMPS_c(&lu->id);
57865e19b50SBarry Smith   if (lu->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1));
579397b6df1SKris Buschelman 
580329ec9b3SHong Zhang   if (lu->size > 1) { /* convert mumps distributed solution to petsc mpi x */
581329ec9b3SHong Zhang     if (!lu->nSolve){ /* create scatter scat_sol */
582329ec9b3SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,lu->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */
583329ec9b3SHong Zhang       for (i=0; i<lu->id.lsol_loc; i++){
584329ec9b3SHong Zhang         lu->id.isol_loc[i] -= 1; /* change Fortran style to C style */
585397b6df1SKris Buschelman       }
58670b3c8c7SBarry Smith       ierr = ISCreateGeneral(PETSC_COMM_SELF,lu->id.lsol_loc,lu->id.isol_loc,PETSC_COPY_VALUES,&is_petsc);CHKERRQ(ierr);  /* to */
587329ec9b3SHong Zhang       ierr = VecScatterCreate(lu->x_seq,is_iden,x,is_petsc,&lu->scat_sol);CHKERRQ(ierr);
5886bf464f9SBarry Smith       ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
5896bf464f9SBarry Smith       ierr = ISDestroy(&is_petsc);CHKERRQ(ierr);
590397b6df1SKris Buschelman     }
591ca9f406cSSatish Balay     ierr = VecScatterBegin(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
592ca9f406cSSatish Balay     ierr = VecScatterEnd(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
593329ec9b3SHong Zhang   }
594329ec9b3SHong Zhang   lu->nSolve++;
595397b6df1SKris Buschelman   PetscFunctionReturn(0);
596397b6df1SKris Buschelman }
597397b6df1SKris Buschelman 
59851d5961aSHong Zhang #undef __FUNCT__
59951d5961aSHong Zhang #define __FUNCT__ "MatSolveTranspose_MUMPS"
60051d5961aSHong Zhang PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x)
60151d5961aSHong Zhang {
60251d5961aSHong Zhang   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
60351d5961aSHong Zhang   PetscErrorCode ierr;
60451d5961aSHong Zhang 
60551d5961aSHong Zhang   PetscFunctionBegin;
60651d5961aSHong Zhang   lu->id.ICNTL(9) = 0;
6070ad0caddSJed Brown   ierr = MatSolve_MUMPS(A,b,x);CHKERRQ(ierr);
60851d5961aSHong Zhang   lu->id.ICNTL(9) = 1;
60951d5961aSHong Zhang   PetscFunctionReturn(0);
61051d5961aSHong Zhang }
61151d5961aSHong Zhang 
612e0b74bf9SHong Zhang #undef __FUNCT__
613e0b74bf9SHong Zhang #define __FUNCT__ "MatMatSolve_MUMPS"
614e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X)
615e0b74bf9SHong Zhang {
616bda8bf91SBarry Smith   PetscErrorCode ierr;
617bda8bf91SBarry Smith   PetscBool      flg;
618bda8bf91SBarry Smith 
619e0b74bf9SHong Zhang   PetscFunctionBegin;
620251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,PETSC_NULL);CHKERRQ(ierr);
621bda8bf91SBarry Smith   if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
622251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,PETSC_NULL);CHKERRQ(ierr);
623bda8bf91SBarry Smith   if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");  SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_MUMPS() is not implemented yet");
624e0b74bf9SHong Zhang   PetscFunctionReturn(0);
625e0b74bf9SHong Zhang }
626e0b74bf9SHong Zhang 
627ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX)
628a58c3f20SHong Zhang /*
629a58c3f20SHong Zhang   input:
630a58c3f20SHong Zhang    F:        numeric factor
631a58c3f20SHong Zhang   output:
632a58c3f20SHong Zhang    nneg:     total number of negative pivots
633a58c3f20SHong Zhang    nzero:    0
634a58c3f20SHong Zhang    npos:     (global dimension of F) - nneg
635a58c3f20SHong Zhang */
636a58c3f20SHong Zhang 
637a58c3f20SHong Zhang #undef __FUNCT__
638a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
639dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
640a58c3f20SHong Zhang {
641a58c3f20SHong Zhang   Mat_MUMPS      *lu =(Mat_MUMPS*)F->spptr;
642dfbe8321SBarry Smith   PetscErrorCode ierr;
643c1490034SHong Zhang   PetscMPIInt    size;
644a58c3f20SHong Zhang 
645a58c3f20SHong Zhang   PetscFunctionBegin;
6467adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)F)->comm,&size);CHKERRQ(ierr);
647bcb30aebSHong 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 */
64865e19b50SBarry Smith   if (size > 1 && lu->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",lu->id.INFOG(13));
649a58c3f20SHong Zhang   if (nneg){
650a58c3f20SHong Zhang     if (!lu->myid){
651a58c3f20SHong Zhang       *nneg = lu->id.INFOG(12);
652a58c3f20SHong Zhang     }
653bcb30aebSHong Zhang     ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr);
654a58c3f20SHong Zhang   }
655a58c3f20SHong Zhang   if (nzero) *nzero = 0;
656d0f46423SBarry Smith   if (npos)  *npos  = F->rmap->N - (*nneg);
657a58c3f20SHong Zhang   PetscFunctionReturn(0);
658a58c3f20SHong Zhang }
659ace3df97SHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */
660a58c3f20SHong Zhang 
661397b6df1SKris Buschelman #undef __FUNCT__
662f6c57405SHong Zhang #define __FUNCT__ "MatFactorNumeric_MUMPS"
6630481f469SBarry Smith PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
664af281ebdSHong Zhang {
665dcd589f8SShri Abhyankar   Mat_MUMPS       *lu =(Mat_MUMPS*)(F)->spptr;
6666849ba73SBarry Smith   PetscErrorCode  ierr;
667e09efc27SHong Zhang   Mat             F_diag;
668ace3abfcSBarry Smith   PetscBool       isMPIAIJ;
669397b6df1SKris Buschelman 
670397b6df1SKris Buschelman   PetscFunctionBegin;
671eb85809bSHong Zhang   ierr = (*lu->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
672397b6df1SKris Buschelman 
673397b6df1SKris Buschelman   /* numerical factorization phase */
674329ec9b3SHong Zhang   /*-------------------------------*/
6753d472b54SHong Zhang   lu->id.job = JOB_FACTNUMERIC;
676958c9bccSBarry Smith   if (!lu->id.ICNTL(18)) {
677a7aca84bSHong Zhang     if (!lu->myid) {
678397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
679*2907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
680*2907cef9SHong Zhang       lu->id.a = (mumps_complex*)lu->val;
681*2907cef9SHong Zhang #else
682397b6df1SKris Buschelman       lu->id.a = (mumps_double_complex*)lu->val;
683*2907cef9SHong Zhang #endif
684397b6df1SKris Buschelman #else
685397b6df1SKris Buschelman       lu->id.a = lu->val;
686397b6df1SKris Buschelman #endif
687397b6df1SKris Buschelman     }
688397b6df1SKris Buschelman   } else {
689397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
690*2907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
691*2907cef9SHong Zhang     lu->id.a_loc = (mumps_complex*)lu->val;
692*2907cef9SHong Zhang #else
693397b6df1SKris Buschelman     lu->id.a_loc = (mumps_double_complex*)lu->val;
694*2907cef9SHong Zhang #endif
695397b6df1SKris Buschelman #else
696397b6df1SKris Buschelman     lu->id.a_loc = lu->val;
697397b6df1SKris Buschelman #endif
698397b6df1SKris Buschelman   }
699*2907cef9SHong Zhang   PetscMUMPS_c(&lu->id);
700397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
70165e19b50SBarry Smith     if (lu->id.INFO(1) == -13) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",lu->id.INFO(2));
70265e19b50SBarry Smith     else SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: INFO(1)=%d, INFO(2)=%d\n",lu->id.INFO(1),lu->id.INFO(2));
703397b6df1SKris Buschelman   }
70465e19b50SBarry Smith   if (!lu->myid && lu->id.ICNTL(16) > 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"  lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16));
705397b6df1SKris Buschelman 
7068ada1bb4SHong Zhang   if (lu->size > 1){
707251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr);
708a214ac2aSShri Abhyankar     if (isMPIAIJ) {
709dcd589f8SShri Abhyankar       F_diag = ((Mat_MPIAIJ *)(F)->data)->A;
710e09efc27SHong Zhang     } else {
711dcd589f8SShri Abhyankar       F_diag = ((Mat_MPISBAIJ *)(F)->data)->A;
712e09efc27SHong Zhang     }
713e09efc27SHong Zhang     F_diag->assembled = PETSC_TRUE;
714329ec9b3SHong Zhang     if (lu->nSolve){
7156bf464f9SBarry Smith       ierr = VecScatterDestroy(&lu->scat_sol);CHKERRQ(ierr);
7160e83c824SBarry Smith       ierr = PetscFree2(lu->id.sol_loc,lu->id.isol_loc);CHKERRQ(ierr);
7176bf464f9SBarry Smith       ierr = VecDestroy(&lu->x_seq);CHKERRQ(ierr);
718329ec9b3SHong Zhang     }
7198ada1bb4SHong Zhang   }
720dcd589f8SShri Abhyankar   (F)->assembled   = PETSC_TRUE;
721397b6df1SKris Buschelman   lu->matstruc     = SAME_NONZERO_PATTERN;
722ace87b0dSHong Zhang   lu->CleanUpMUMPS = PETSC_TRUE;
723329ec9b3SHong Zhang   lu->nSolve       = 0;
72467877ebaSShri Abhyankar 
72567877ebaSShri Abhyankar   if (lu->size > 1){
72667877ebaSShri Abhyankar     /* distributed solution */
72767877ebaSShri Abhyankar     if (!lu->nSolve){
72867877ebaSShri Abhyankar       /* Create x_seq=sol_loc for repeated use */
72967877ebaSShri Abhyankar       PetscInt    lsol_loc;
73067877ebaSShri Abhyankar       PetscScalar *sol_loc;
73167877ebaSShri Abhyankar       lsol_loc = lu->id.INFO(23); /* length of sol_loc */
73267877ebaSShri Abhyankar       ierr = PetscMalloc2(lsol_loc,PetscScalar,&sol_loc,lsol_loc,PetscInt,&lu->id.isol_loc);CHKERRQ(ierr);
73367877ebaSShri Abhyankar       lu->id.lsol_loc = lsol_loc;
73467877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
735*2907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
736*2907cef9SHong Zhang       lu->id.sol_loc  = (mumps_complex*)sol_loc;
737*2907cef9SHong Zhang #else
73867877ebaSShri Abhyankar       lu->id.sol_loc  = (mumps_double_complex*)sol_loc;
739*2907cef9SHong Zhang #endif
74067877ebaSShri Abhyankar #else
74167877ebaSShri Abhyankar       lu->id.sol_loc  = sol_loc;
74267877ebaSShri Abhyankar #endif
743778a2246SBarry Smith       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&lu->x_seq);CHKERRQ(ierr);
74467877ebaSShri Abhyankar     }
74567877ebaSShri Abhyankar   }
746397b6df1SKris Buschelman   PetscFunctionReturn(0);
747397b6df1SKris Buschelman }
748397b6df1SKris Buschelman 
7499a2535b5SHong Zhang /* Sets MUMPS options from the options database */
750dcd589f8SShri Abhyankar #undef __FUNCT__
7519a2535b5SHong Zhang #define __FUNCT__ "PetscSetMUMPSFromOptions"
7529a2535b5SHong Zhang PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A)
753dcd589f8SShri Abhyankar {
7549a2535b5SHong Zhang   Mat_MUMPS        *mumps = (Mat_MUMPS*)F->spptr;
755dcd589f8SShri Abhyankar   PetscErrorCode   ierr;
756dcd589f8SShri Abhyankar   PetscInt         icntl;
757ace3abfcSBarry Smith   PetscBool        flg;
758dcd589f8SShri Abhyankar 
759dcd589f8SShri Abhyankar   PetscFunctionBegin;
760dcd589f8SShri Abhyankar   ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
7619a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);CHKERRQ(ierr);
7629a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(1) = icntl;
7639a2535b5SHong 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);
7649a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(2) = icntl;
7659a2535b5SHong 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);
7669a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(3) = icntl;
767dcd589f8SShri Abhyankar 
7689a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",mumps->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
7699a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(4) = icntl;
7709a2535b5SHong Zhang   if (mumps->id.ICNTL(4) || PetscLogPrintInfo ) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */
7719a2535b5SHong Zhang 
7729a2535b5SHong 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);
7739a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(6) = icntl;
7749a2535b5SHong Zhang 
7759a2535b5SHong 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);
776dcd589f8SShri Abhyankar   if (flg) {
7779a2535b5SHong Zhang     if (icntl== 1 && mumps->size > 1){
778e32f2f54SBarry Smith       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");
779dcd589f8SShri Abhyankar     } else {
7809a2535b5SHong Zhang       mumps->id.ICNTL(7) = icntl;
781dcd589f8SShri Abhyankar     }
782dcd589f8SShri Abhyankar   }
783e0b74bf9SHong Zhang 
78470544d5fSHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_8","ICNTL(8): scaling strategy (-2 to 8 or 77)","None",mumps->id.ICNTL(8),&mumps->id.ICNTL(8),PETSC_NULL);CHKERRQ(ierr);
7859a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",mumps->id.ICNTL(10),&mumps->id.ICNTL(10),PETSC_NULL);CHKERRQ(ierr);
7869a2535b5SHong Zhang   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),PETSC_NULL);CHKERRQ(ierr);
78770544d5fSHong Zhang   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),PETSC_NULL);CHKERRQ(ierr);
7889a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control: with or without ScaLAPACK","None",mumps->id.ICNTL(13),&mumps->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr);
7899a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage of estimated workspace increase","None",mumps->id.ICNTL(14),&mumps->id.ICNTL(14),PETSC_NULL);CHKERRQ(ierr);
7909a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): Schur complement","None",mumps->id.ICNTL(19),&mumps->id.ICNTL(19),PETSC_NULL);CHKERRQ(ierr);
7919a2535b5SHong Zhang 
7929a2535b5SHong Zhang   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),PETSC_NULL);CHKERRQ(ierr);
7939a2535b5SHong Zhang   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),PETSC_NULL);CHKERRQ(ierr);
7949a2535b5SHong Zhang   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),PETSC_NULL);CHKERRQ(ierr);
7959a2535b5SHong Zhang   if (mumps->id.ICNTL(24)){
7969a2535b5SHong Zhang     mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */
797d7ebd59bSHong Zhang   }
798d7ebd59bSHong Zhang 
7999a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_25","ICNTL(25): computation of a null space basis","None",mumps->id.ICNTL(25),&mumps->id.ICNTL(25),PETSC_NULL);CHKERRQ(ierr);
8009a2535b5SHong Zhang   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),PETSC_NULL);CHKERRQ(ierr);
8019a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_27","ICNTL(27): experimental parameter","None",mumps->id.ICNTL(27),&mumps->id.ICNTL(27),PETSC_NULL);CHKERRQ(ierr);
8029a2535b5SHong Zhang   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),PETSC_NULL);CHKERRQ(ierr);
8039a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_29","ICNTL(29): parallel ordering 1 = ptscotch 2 = parmetis","None",mumps->id.ICNTL(29),&mumps->id.ICNTL(29),PETSC_NULL);CHKERRQ(ierr);
8049a2535b5SHong Zhang   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),PETSC_NULL);CHKERRQ(ierr);
8059a2535b5SHong Zhang   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),PETSC_NULL);CHKERRQ(ierr);
8069a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),PETSC_NULL);CHKERRQ(ierr);
807dcd589f8SShri Abhyankar 
8089a2535b5SHong Zhang   ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),PETSC_NULL);CHKERRQ(ierr);
8099a2535b5SHong Zhang   ierr = PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",mumps->id.CNTL(2),&mumps->id.CNTL(2),PETSC_NULL);CHKERRQ(ierr);
8109a2535b5SHong Zhang   ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),PETSC_NULL);CHKERRQ(ierr);
8119a2535b5SHong Zhang   ierr = PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",mumps->id.CNTL(4),&mumps->id.CNTL(4),PETSC_NULL);CHKERRQ(ierr);
8129a2535b5SHong Zhang   ierr = PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",mumps->id.CNTL(5),&mumps->id.CNTL(5),PETSC_NULL);CHKERRQ(ierr);
813e5bb22a1SHong Zhang 
814e5bb22a1SHong Zhang   ierr = PetscOptionsString("-mat_mumps_ooc_tmpdir", "out of core directory", "None", mumps->id.ooc_tmpdir, mumps->id.ooc_tmpdir, 256, PETSC_NULL);
815dcd589f8SShri Abhyankar   PetscOptionsEnd();
816dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
817dcd589f8SShri Abhyankar }
818dcd589f8SShri Abhyankar 
819dcd589f8SShri Abhyankar #undef __FUNCT__
820dcd589f8SShri Abhyankar #define __FUNCT__ "PetscInitializeMUMPS"
821f697e70eSHong Zhang PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS* mumps)
822dcd589f8SShri Abhyankar {
823dcd589f8SShri Abhyankar   PetscErrorCode  ierr;
824dcd589f8SShri Abhyankar 
825dcd589f8SShri Abhyankar   PetscFunctionBegin;
826f697e70eSHong Zhang   ierr = MPI_Comm_rank(((PetscObject)A)->comm, &mumps->myid);
827f697e70eSHong Zhang   ierr = MPI_Comm_size(((PetscObject)A)->comm,&mumps->size);CHKERRQ(ierr);
828f697e70eSHong Zhang   ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(mumps->comm_mumps));CHKERRQ(ierr);
829f697e70eSHong Zhang   mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps);
830f697e70eSHong Zhang 
831f697e70eSHong Zhang   mumps->id.job = JOB_INIT;
832f697e70eSHong Zhang   mumps->id.par = 1;  /* host participates factorizaton and solve */
833f697e70eSHong Zhang   mumps->id.sym = mumps->sym;
834*2907cef9SHong Zhang   PetscMUMPS_c(&mumps->id);
835f697e70eSHong Zhang 
836f697e70eSHong Zhang   mumps->CleanUpMUMPS = PETSC_FALSE;
837f697e70eSHong Zhang   mumps->scat_rhs     = PETSC_NULL;
838f697e70eSHong Zhang   mumps->scat_sol     = PETSC_NULL;
839f697e70eSHong Zhang   mumps->nSolve       = 0;
8409a2535b5SHong Zhang 
84170544d5fSHong Zhang   /* set PETSc-MUMPS default options - override MUMPS default */
8429a2535b5SHong Zhang   mumps->id.ICNTL(3) = 0;
8439a2535b5SHong Zhang   mumps->id.ICNTL(4) = 0;
8449a2535b5SHong Zhang   if (mumps->size == 1){
8459a2535b5SHong Zhang     mumps->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
8469a2535b5SHong Zhang   } else {
8479a2535b5SHong Zhang     mumps->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
84870544d5fSHong Zhang     mumps->id.ICNTL(21) = 1;   /* distributed solution */
8499a2535b5SHong Zhang   }
850dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
851dcd589f8SShri Abhyankar }
852dcd589f8SShri Abhyankar 
853397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */
854397b6df1SKris Buschelman #undef __FUNCT__
855f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
8560481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
857b24902e0SBarry Smith {
858719d5645SBarry Smith   Mat_MUMPS          *lu = (Mat_MUMPS*)F->spptr;
859dcd589f8SShri Abhyankar   PetscErrorCode     ierr;
86067877ebaSShri Abhyankar   Vec                b;
86167877ebaSShri Abhyankar   IS                 is_iden;
86267877ebaSShri Abhyankar   const PetscInt     M = A->rmap->N;
863397b6df1SKris Buschelman 
864397b6df1SKris Buschelman   PetscFunctionBegin;
865b24902e0SBarry Smith   lu->matstruc = DIFFERENT_NONZERO_PATTERN;
866dcd589f8SShri Abhyankar 
8679a2535b5SHong Zhang   /* Set MUMPS options from the options database */
8689a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
869dcd589f8SShri Abhyankar 
870eb85809bSHong Zhang   ierr = (*lu->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
871dcd589f8SShri Abhyankar 
87267877ebaSShri Abhyankar   /* analysis phase */
87367877ebaSShri Abhyankar   /*----------------*/
8743d472b54SHong Zhang   lu->id.job = JOB_FACTSYMBOLIC;
87567877ebaSShri Abhyankar   lu->id.n = M;
87667877ebaSShri Abhyankar   switch (lu->id.ICNTL(18)){
87767877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
87867877ebaSShri Abhyankar     if (!lu->myid) {
87967877ebaSShri Abhyankar       lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
88067877ebaSShri Abhyankar       if (lu->id.ICNTL(6)>1){
88167877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
882*2907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
883*2907cef9SHong Zhang         lu->id.a = (mumps_complex*)lu->val;
884*2907cef9SHong Zhang #else
88567877ebaSShri Abhyankar         lu->id.a = (mumps_double_complex*)lu->val;
886*2907cef9SHong Zhang #endif
88767877ebaSShri Abhyankar #else
88867877ebaSShri Abhyankar         lu->id.a = lu->val;
88967877ebaSShri Abhyankar #endif
89067877ebaSShri Abhyankar       }
891e0b74bf9SHong Zhang       if (lu->id.ICNTL(7) == 1){ /* use user-provide matrix ordering */
892e0b74bf9SHong Zhang         if (!lu->myid) {
893e0b74bf9SHong Zhang           const PetscInt *idx;
894e0b74bf9SHong Zhang           PetscInt i,*perm_in;
895e0b74bf9SHong Zhang           ierr = PetscMalloc(M*sizeof(PetscInt),&perm_in);CHKERRQ(ierr);
896e0b74bf9SHong Zhang           ierr = ISGetIndices(r,&idx);CHKERRQ(ierr);
897e0b74bf9SHong Zhang           lu->id.perm_in = perm_in;
898e0b74bf9SHong Zhang           for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */
899e0b74bf9SHong Zhang           ierr = ISRestoreIndices(r,&idx);CHKERRQ(ierr);
900e0b74bf9SHong Zhang         }
901e0b74bf9SHong Zhang       }
90267877ebaSShri Abhyankar     }
90367877ebaSShri Abhyankar     break;
90467877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
90567877ebaSShri Abhyankar     lu->id.nz_loc = lu->nz;
90667877ebaSShri Abhyankar     lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
90767877ebaSShri Abhyankar     if (lu->id.ICNTL(6)>1) {
90867877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
909*2907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
910*2907cef9SHong Zhang       lu->id.a_loc = (mumps_complex*)lu->val;
911*2907cef9SHong Zhang #else
91267877ebaSShri Abhyankar       lu->id.a_loc = (mumps_double_complex*)lu->val;
913*2907cef9SHong Zhang #endif
91467877ebaSShri Abhyankar #else
91567877ebaSShri Abhyankar       lu->id.a_loc = lu->val;
91667877ebaSShri Abhyankar #endif
91767877ebaSShri Abhyankar     }
91867877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
91967877ebaSShri Abhyankar     if (!lu->myid){
92067877ebaSShri Abhyankar       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr);
92167877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
92267877ebaSShri Abhyankar     } else {
92367877ebaSShri Abhyankar       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr);
92467877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
92567877ebaSShri Abhyankar     }
92667877ebaSShri Abhyankar     ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr);
92767877ebaSShri Abhyankar     ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
92867877ebaSShri Abhyankar     ierr = VecSetFromOptions(b);CHKERRQ(ierr);
92967877ebaSShri Abhyankar 
93067877ebaSShri Abhyankar     ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr);
9316bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
9326bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
93367877ebaSShri Abhyankar     break;
93467877ebaSShri Abhyankar     }
935*2907cef9SHong Zhang   PetscMUMPS_c(&lu->id);
93667877ebaSShri Abhyankar   if (lu->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
93767877ebaSShri Abhyankar 
938719d5645SBarry Smith   F->ops->lufactornumeric  = MatFactorNumeric_MUMPS;
939dcd589f8SShri Abhyankar   F->ops->solve            = MatSolve_MUMPS;
94051d5961aSHong Zhang   F->ops->solvetranspose   = MatSolveTranspose_MUMPS;
941e0b74bf9SHong Zhang   F->ops->matsolve         = MatMatSolve_MUMPS;
942b24902e0SBarry Smith   PetscFunctionReturn(0);
943b24902e0SBarry Smith }
944b24902e0SBarry Smith 
945450b117fSShri Abhyankar /* Note the Petsc r and c permutations are ignored */
946450b117fSShri Abhyankar #undef __FUNCT__
947450b117fSShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_BAIJMUMPS"
948450b117fSShri Abhyankar PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
949450b117fSShri Abhyankar {
950dcd589f8SShri Abhyankar 
951450b117fSShri Abhyankar   Mat_MUMPS       *lu = (Mat_MUMPS*)F->spptr;
952dcd589f8SShri Abhyankar   PetscErrorCode  ierr;
95367877ebaSShri Abhyankar   Vec             b;
95467877ebaSShri Abhyankar   IS              is_iden;
95567877ebaSShri Abhyankar   const PetscInt  M = A->rmap->N;
956450b117fSShri Abhyankar 
957450b117fSShri Abhyankar   PetscFunctionBegin;
958450b117fSShri Abhyankar   lu->matstruc = DIFFERENT_NONZERO_PATTERN;
959dcd589f8SShri Abhyankar 
9609a2535b5SHong Zhang   /* Set MUMPS options from the options database */
9619a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
962dcd589f8SShri Abhyankar 
963eb85809bSHong Zhang   ierr = (*lu->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
96467877ebaSShri Abhyankar 
96567877ebaSShri Abhyankar   /* analysis phase */
96667877ebaSShri Abhyankar   /*----------------*/
9673d472b54SHong Zhang   lu->id.job = JOB_FACTSYMBOLIC;
96867877ebaSShri Abhyankar   lu->id.n = M;
96967877ebaSShri Abhyankar   switch (lu->id.ICNTL(18)){
97067877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
97167877ebaSShri Abhyankar     if (!lu->myid) {
97267877ebaSShri Abhyankar       lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
97367877ebaSShri Abhyankar       if (lu->id.ICNTL(6)>1){
97467877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
975*2907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
976*2907cef9SHong Zhang         lu->id.a = (mumps_complex*)lu->val;
977*2907cef9SHong Zhang #else
97867877ebaSShri Abhyankar         lu->id.a = (mumps_double_complex*)lu->val;
979*2907cef9SHong Zhang #endif
98067877ebaSShri Abhyankar #else
98167877ebaSShri Abhyankar         lu->id.a = lu->val;
98267877ebaSShri Abhyankar #endif
98367877ebaSShri Abhyankar       }
98467877ebaSShri Abhyankar     }
98567877ebaSShri Abhyankar     break;
98667877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
98767877ebaSShri Abhyankar     lu->id.nz_loc = lu->nz;
98867877ebaSShri Abhyankar     lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
98967877ebaSShri Abhyankar     if (lu->id.ICNTL(6)>1) {
99067877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
991*2907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
992*2907cef9SHong Zhang       lu->id.a_loc = (mumps_complex*)lu->val;
993*2907cef9SHong Zhang #else
99467877ebaSShri Abhyankar       lu->id.a_loc = (mumps_double_complex*)lu->val;
995*2907cef9SHong Zhang #endif
99667877ebaSShri Abhyankar #else
99767877ebaSShri Abhyankar       lu->id.a_loc = lu->val;
99867877ebaSShri Abhyankar #endif
99967877ebaSShri Abhyankar     }
100067877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
100167877ebaSShri Abhyankar     if (!lu->myid){
100267877ebaSShri Abhyankar       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr);
100367877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
100467877ebaSShri Abhyankar     } else {
100567877ebaSShri Abhyankar       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr);
100667877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
100767877ebaSShri Abhyankar     }
100867877ebaSShri Abhyankar     ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr);
100967877ebaSShri Abhyankar     ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
101067877ebaSShri Abhyankar     ierr = VecSetFromOptions(b);CHKERRQ(ierr);
101167877ebaSShri Abhyankar 
101267877ebaSShri Abhyankar     ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr);
10136bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
10146bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
101567877ebaSShri Abhyankar     break;
101667877ebaSShri Abhyankar     }
1017*2907cef9SHong Zhang   PetscMUMPS_c(&lu->id);
101867877ebaSShri Abhyankar   if (lu->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
101967877ebaSShri Abhyankar 
1020450b117fSShri Abhyankar   F->ops->lufactornumeric  = MatFactorNumeric_MUMPS;
1021dcd589f8SShri Abhyankar   F->ops->solve            = MatSolve_MUMPS;
102251d5961aSHong Zhang   F->ops->solvetranspose   = MatSolveTranspose_MUMPS;
1023450b117fSShri Abhyankar   PetscFunctionReturn(0);
1024450b117fSShri Abhyankar }
1025b24902e0SBarry Smith 
1026141f4205SHong Zhang /* Note the Petsc r permutation and factor info are ignored */
1027397b6df1SKris Buschelman #undef __FUNCT__
102867877ebaSShri Abhyankar #define __FUNCT__ "MatCholeskyFactorSymbolic_MUMPS"
102967877ebaSShri Abhyankar PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
1030b24902e0SBarry Smith {
10312792810eSHong Zhang   Mat_MUMPS          *lu = (Mat_MUMPS*)F->spptr;
1032dcd589f8SShri Abhyankar   PetscErrorCode     ierr;
103367877ebaSShri Abhyankar   Vec                b;
103467877ebaSShri Abhyankar   IS                 is_iden;
103567877ebaSShri Abhyankar   const PetscInt     M = A->rmap->N;
1036397b6df1SKris Buschelman 
1037397b6df1SKris Buschelman   PetscFunctionBegin;
1038b24902e0SBarry Smith   lu->matstruc = DIFFERENT_NONZERO_PATTERN;
1039dcd589f8SShri Abhyankar 
10409a2535b5SHong Zhang   /* Set MUMPS options from the options database */
10419a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1042dcd589f8SShri Abhyankar 
1043eb85809bSHong Zhang   ierr = (*lu->ConvertToTriples)(A, 1 , MAT_INITIAL_MATRIX, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
1044dcd589f8SShri Abhyankar 
104567877ebaSShri Abhyankar   /* analysis phase */
104667877ebaSShri Abhyankar   /*----------------*/
10473d472b54SHong Zhang   lu->id.job = JOB_FACTSYMBOLIC;
104867877ebaSShri Abhyankar   lu->id.n = M;
104967877ebaSShri Abhyankar   switch (lu->id.ICNTL(18)){
105067877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
105167877ebaSShri Abhyankar     if (!lu->myid) {
105267877ebaSShri Abhyankar       lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
105367877ebaSShri Abhyankar       if (lu->id.ICNTL(6)>1){
105467877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
1055*2907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1056*2907cef9SHong Zhang         lu->id.a = (mumps_complex*)lu->val;
1057*2907cef9SHong Zhang #else
105867877ebaSShri Abhyankar         lu->id.a = (mumps_double_complex*)lu->val;
1059*2907cef9SHong Zhang #endif
106067877ebaSShri Abhyankar #else
106167877ebaSShri Abhyankar         lu->id.a = lu->val;
106267877ebaSShri Abhyankar #endif
106367877ebaSShri Abhyankar       }
106467877ebaSShri Abhyankar     }
106567877ebaSShri Abhyankar     break;
106667877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
106767877ebaSShri Abhyankar     lu->id.nz_loc = lu->nz;
106867877ebaSShri Abhyankar     lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
106967877ebaSShri Abhyankar     if (lu->id.ICNTL(6)>1) {
107067877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
1071*2907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1072*2907cef9SHong Zhang       lu->id.a_loc = (mumps_complex*)lu->val;
1073*2907cef9SHong Zhang #else
107467877ebaSShri Abhyankar       lu->id.a_loc = (mumps_double_complex*)lu->val;
1075*2907cef9SHong Zhang #endif
107667877ebaSShri Abhyankar #else
107767877ebaSShri Abhyankar       lu->id.a_loc = lu->val;
107867877ebaSShri Abhyankar #endif
107967877ebaSShri Abhyankar     }
108067877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
108167877ebaSShri Abhyankar     if (!lu->myid){
108267877ebaSShri Abhyankar       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr);
108367877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
108467877ebaSShri Abhyankar     } else {
108567877ebaSShri Abhyankar       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr);
108667877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
108767877ebaSShri Abhyankar     }
108867877ebaSShri Abhyankar     ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr);
108967877ebaSShri Abhyankar     ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
109067877ebaSShri Abhyankar     ierr = VecSetFromOptions(b);CHKERRQ(ierr);
109167877ebaSShri Abhyankar 
109267877ebaSShri Abhyankar     ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr);
10936bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
10946bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
109567877ebaSShri Abhyankar     break;
109667877ebaSShri Abhyankar     }
1097*2907cef9SHong Zhang   PetscMUMPS_c(&lu->id);
109867877ebaSShri Abhyankar   if (lu->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
109967877ebaSShri Abhyankar 
11002792810eSHong Zhang   F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
1101dcd589f8SShri Abhyankar   F->ops->solve                 = MatSolve_MUMPS;
110251d5961aSHong Zhang   F->ops->solvetranspose        = MatSolve_MUMPS;
1103377a2891SHong Zhang   F->ops->matsolve              = MatMatSolve_MUMPS;
1104db4efbfdSBarry Smith #if !defined(PETSC_USE_COMPLEX)
110505aa0992SJose Roman   F->ops->getinertia            = MatGetInertia_SBAIJMUMPS;
110605aa0992SJose Roman #else
110705aa0992SJose Roman   F->ops->getinertia            = PETSC_NULL;
1108db4efbfdSBarry Smith #endif
1109b24902e0SBarry Smith   PetscFunctionReturn(0);
1110b24902e0SBarry Smith }
1111b24902e0SBarry Smith 
1112397b6df1SKris Buschelman #undef __FUNCT__
111364e6c443SBarry Smith #define __FUNCT__ "MatView_MUMPS"
111464e6c443SBarry Smith PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
111574ed9c26SBarry Smith {
1116f6c57405SHong Zhang   PetscErrorCode    ierr;
111764e6c443SBarry Smith   PetscBool         iascii;
111864e6c443SBarry Smith   PetscViewerFormat format;
111964e6c443SBarry Smith   Mat_MUMPS         *lu=(Mat_MUMPS*)A->spptr;
1120f6c57405SHong Zhang 
1121f6c57405SHong Zhang   PetscFunctionBegin;
112264e6c443SBarry Smith   /* check if matrix is mumps type */
112364e6c443SBarry Smith   if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0);
112464e6c443SBarry Smith 
1125251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
112664e6c443SBarry Smith   if (iascii) {
112764e6c443SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
112864e6c443SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO){
112964e6c443SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
113064e6c443SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                   %d \n",lu->id.sym);CHKERRQ(ierr);
113164e6c443SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):            %d \n",lu->id.par);CHKERRQ(ierr);
1132f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):         %d \n",lu->id.ICNTL(1));CHKERRQ(ierr);
1133f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg): %d \n",lu->id.ICNTL(2));CHKERRQ(ierr);
1134f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):   %d \n",lu->id.ICNTL(3));CHKERRQ(ierr);
1135f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):        %d \n",lu->id.ICNTL(4));CHKERRQ(ierr);
1136f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):         %d \n",lu->id.ICNTL(5));CHKERRQ(ierr);
1137f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):        %d \n",lu->id.ICNTL(6));CHKERRQ(ierr);
1138d06efebcSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (sequentia matrix ordering):%d \n",lu->id.ICNTL(7));CHKERRQ(ierr);
1139f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scalling strategy):        %d \n",lu->id.ICNTL(8));CHKERRQ(ierr);
1140f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements):  %d \n",lu->id.ICNTL(10));CHKERRQ(ierr);
1141f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):          %d \n",lu->id.ICNTL(11));CHKERRQ(ierr);
114234ed7027SBarry Smith       if (lu->id.ICNTL(11)>0) {
114334ed7027SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(4) (inf norm of input mat):        %g\n",lu->id.RINFOG(4));CHKERRQ(ierr);
114434ed7027SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(5) (inf norm of solution):         %g\n",lu->id.RINFOG(5));CHKERRQ(ierr);
114534ed7027SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(6) (inf norm of residual):         %g\n",lu->id.RINFOG(6));CHKERRQ(ierr);
114634ed7027SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",lu->id.RINFOG(7),lu->id.RINFOG(8));CHKERRQ(ierr);
114734ed7027SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(9) (error estimate):               %g \n",lu->id.RINFOG(9));CHKERRQ(ierr);
114834ed7027SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",lu->id.RINFOG(10),lu->id.RINFOG(11));CHKERRQ(ierr);
1149f6c57405SHong Zhang       }
1150f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",lu->id.ICNTL(12));CHKERRQ(ierr);
1151f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",lu->id.ICNTL(13));CHKERRQ(ierr);
1152f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr);
1153f6c57405SHong Zhang       /* ICNTL(15-17) not used */
1154f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",lu->id.ICNTL(18));CHKERRQ(ierr);
1155f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Shur complement info):                       %d \n",lu->id.ICNTL(19));CHKERRQ(ierr);
1156f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",lu->id.ICNTL(20));CHKERRQ(ierr);
1157f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d \n",lu->id.ICNTL(21));CHKERRQ(ierr);
1158c0165424SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",lu->id.ICNTL(22));CHKERRQ(ierr);
1159c0165424SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",lu->id.ICNTL(23));CHKERRQ(ierr);
1160c0165424SHong Zhang 
1161c0165424SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",lu->id.ICNTL(24));CHKERRQ(ierr);
1162c0165424SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",lu->id.ICNTL(25));CHKERRQ(ierr);
1163c0165424SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",lu->id.ICNTL(26));CHKERRQ(ierr);
1164c0165424SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",lu->id.ICNTL(27));CHKERRQ(ierr);
116542179a6aSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(28) (use parallel or sequential ordering):        %d \n",lu->id.ICNTL(28));CHKERRQ(ierr);
116642179a6aSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(29) (parallel ordering):                          %d \n",lu->id.ICNTL(29));CHKERRQ(ierr);
116742179a6aSHong Zhang 
116842179a6aSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(30) (user-specified set of entries in inv(A)):    %d \n",lu->id.ICNTL(30));CHKERRQ(ierr);
116942179a6aSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(31) (factors is discarded in the solve phase):    %d \n",lu->id.ICNTL(31));CHKERRQ(ierr);
117042179a6aSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(33) (compute determinant):                        %d \n",lu->id.ICNTL(33));CHKERRQ(ierr);
1171f6c57405SHong Zhang 
1172f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",lu->id.CNTL(1));CHKERRQ(ierr);
1173f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr);
1174f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",lu->id.CNTL(3));CHKERRQ(ierr);
1175f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",lu->id.CNTL(4));CHKERRQ(ierr);
1176c0165424SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",lu->id.CNTL(5));CHKERRQ(ierr);
1177f6c57405SHong Zhang 
1178f6c57405SHong Zhang       /* infomation local to each processor */
117934ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);
11807b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
118134ed7027SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr);
118234ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
118334ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);
118434ed7027SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr);
118534ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
118634ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);
118734ed7027SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr);
118834ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1189f6c57405SHong Zhang 
119034ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr);
119134ed7027SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] %d \n",lu->myid,lu->id.INFO(15));CHKERRQ(ierr);
119234ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1193f6c57405SHong Zhang 
119434ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr);
119534ed7027SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",lu->myid,lu->id.INFO(16));CHKERRQ(ierr);
119634ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1197f6c57405SHong Zhang 
119834ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);
119934ed7027SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",lu->myid,lu->id.INFO(23));CHKERRQ(ierr);
120034ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
12017b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
1202f6c57405SHong Zhang 
1203f6c57405SHong Zhang       if (!lu->myid){ /* information from the host */
1204f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr);
1205f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr);
1206f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr);
120742179a6aSHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  (RINFOG(12) RINFOG(13))*2^INFOG(34) (determinant): (%g,%g)*(2^%d)\n",lu->id.RINFOG(12),lu->id.RINFOG(13),lu->id.INFOG(34));CHKERRQ(ierr);
1208f6c57405SHong Zhang 
1209f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr);
1210f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr);
1211f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr);
1212f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr);
12132bd8dccdSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively use after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr);
1214f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr);
1215f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",lu->id.INFOG(9));CHKERRQ(ierr);
1216f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",lu->id.INFOG(10));CHKERRQ(ierr);
1217f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",lu->id.INFOG(11));CHKERRQ(ierr);
1218f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr);
1219f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr);
1220f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr);
1221f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr);
1222f6c57405SHong 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",lu->id.INFOG(16));CHKERRQ(ierr);
1223f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(17) (estimated size of all MUMPS internal data for factorization after analysis: sum over all processors): %d \n",lu->id.INFOG(17));CHKERRQ(ierr);
1224f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(18) (size of all MUMPS internal data allocated during factorization: value on the most memory consuming processor): %d \n",lu->id.INFOG(18));CHKERRQ(ierr);
1225f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(19) (size of all MUMPS internal data allocated during factorization: sum over all processors): %d \n",lu->id.INFOG(19));CHKERRQ(ierr);
1226f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr);
1227f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(21) (size in MB of memory effectively used during factorization - value on the most memory consuming processor): %d \n",lu->id.INFOG(21));CHKERRQ(ierr);
1228f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(22) (size in MB of memory effectively used during factorization - sum over all processors): %d \n",lu->id.INFOG(22));CHKERRQ(ierr);
1229f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",lu->id.INFOG(23));CHKERRQ(ierr);
1230f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",lu->id.INFOG(24));CHKERRQ(ierr);
1231f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",lu->id.INFOG(25));CHKERRQ(ierr);
1232f6c57405SHong Zhang       }
1233f6c57405SHong Zhang     }
1234cb828f0fSHong Zhang   }
1235f6c57405SHong Zhang   PetscFunctionReturn(0);
1236f6c57405SHong Zhang }
1237f6c57405SHong Zhang 
123835bd34faSBarry Smith #undef __FUNCT__
123935bd34faSBarry Smith #define __FUNCT__ "MatGetInfo_MUMPS"
124035bd34faSBarry Smith PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
124135bd34faSBarry Smith {
1242cb828f0fSHong Zhang   Mat_MUMPS  *mumps =(Mat_MUMPS*)A->spptr;
124335bd34faSBarry Smith 
124435bd34faSBarry Smith   PetscFunctionBegin;
124535bd34faSBarry Smith   info->block_size        = 1.0;
1246cb828f0fSHong Zhang   info->nz_allocated      = mumps->id.INFOG(20);
1247cb828f0fSHong Zhang   info->nz_used           = mumps->id.INFOG(20);
124835bd34faSBarry Smith   info->nz_unneeded       = 0.0;
124935bd34faSBarry Smith   info->assemblies        = 0.0;
125035bd34faSBarry Smith   info->mallocs           = 0.0;
125135bd34faSBarry Smith   info->memory            = 0.0;
125235bd34faSBarry Smith   info->fill_ratio_given  = 0;
125335bd34faSBarry Smith   info->fill_ratio_needed = 0;
125435bd34faSBarry Smith   info->factor_mallocs    = 0;
125535bd34faSBarry Smith   PetscFunctionReturn(0);
125635bd34faSBarry Smith }
125735bd34faSBarry Smith 
12585ccb76cbSHong Zhang /* -------------------------------------------------------------------------------------------*/
12595ccb76cbSHong Zhang #undef __FUNCT__
12605ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl_MUMPS"
12615ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival)
12625ccb76cbSHong Zhang {
12635ccb76cbSHong Zhang   Mat_MUMPS *lu =(Mat_MUMPS*)F->spptr;
12645ccb76cbSHong Zhang 
12655ccb76cbSHong Zhang   PetscFunctionBegin;
12665ccb76cbSHong Zhang   lu->id.ICNTL(icntl) = ival;
12675ccb76cbSHong Zhang   PetscFunctionReturn(0);
12685ccb76cbSHong Zhang }
12695ccb76cbSHong Zhang 
12705ccb76cbSHong Zhang #undef __FUNCT__
12715ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl"
12725ccb76cbSHong Zhang /*@
12735ccb76cbSHong Zhang   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()
12745ccb76cbSHong Zhang 
12755ccb76cbSHong Zhang    Logically Collective on Mat
12765ccb76cbSHong Zhang 
12775ccb76cbSHong Zhang    Input Parameters:
12785ccb76cbSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
12795ccb76cbSHong Zhang .  icntl - index of MUMPS parameter array ICNTL()
12805ccb76cbSHong Zhang -  ival - value of MUMPS ICNTL(icntl)
12815ccb76cbSHong Zhang 
12825ccb76cbSHong Zhang   Options Database:
12835ccb76cbSHong Zhang .   -mat_mumps_icntl_<icntl> <ival>
12845ccb76cbSHong Zhang 
12855ccb76cbSHong Zhang    Level: beginner
12865ccb76cbSHong Zhang 
12875ccb76cbSHong Zhang    References: MUMPS Users' Guide
12885ccb76cbSHong Zhang 
12895ccb76cbSHong Zhang .seealso: MatGetFactor()
12905ccb76cbSHong Zhang @*/
12915ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival)
12925ccb76cbSHong Zhang {
12935ccb76cbSHong Zhang   PetscErrorCode ierr;
12945ccb76cbSHong Zhang 
12955ccb76cbSHong Zhang   PetscFunctionBegin;
12965ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
12975ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,ival,3);
12985ccb76cbSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr);
12995ccb76cbSHong Zhang   PetscFunctionReturn(0);
13005ccb76cbSHong Zhang }
13015ccb76cbSHong Zhang 
130224b6179bSKris Buschelman /*MC
13032692d6eeSBarry Smith   MATSOLVERMUMPS -  A matrix type providing direct solvers (LU and Cholesky) for
130424b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
130524b6179bSKris Buschelman 
130641c8de11SBarry Smith   Works with MATAIJ and MATSBAIJ matrices
130724b6179bSKris Buschelman 
130824b6179bSKris Buschelman   Options Database Keys:
1309fb8376fbSHong Zhang + -mat_mumps_icntl_4 <0,...,4> - print level
131024b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
131164e6c443SBarry Smith . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guidec)
131224b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
131324b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
131494b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
131524b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
131624b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
131724b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
131824b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
131924b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
132024b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
132124b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
132224b6179bSKris Buschelman 
132324b6179bSKris Buschelman   Level: beginner
132424b6179bSKris Buschelman 
132541c8de11SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage
132641c8de11SBarry Smith 
132724b6179bSKris Buschelman M*/
132824b6179bSKris Buschelman 
13292877fffaSHong Zhang EXTERN_C_BEGIN
133035bd34faSBarry Smith #undef __FUNCT__
133135bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_mumps"
133235bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type)
133335bd34faSBarry Smith {
133435bd34faSBarry Smith   PetscFunctionBegin;
13352692d6eeSBarry Smith   *type = MATSOLVERMUMPS;
133635bd34faSBarry Smith   PetscFunctionReturn(0);
133735bd34faSBarry Smith }
133835bd34faSBarry Smith EXTERN_C_END
133935bd34faSBarry Smith 
134035bd34faSBarry Smith EXTERN_C_BEGIN
1341bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI AIJ matrices */
13422877fffaSHong Zhang #undef __FUNCT__
1343bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_aij_mumps"
1344bccb9932SShri Abhyankar PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
13452877fffaSHong Zhang {
13462877fffaSHong Zhang   Mat            B;
13472877fffaSHong Zhang   PetscErrorCode ierr;
13482877fffaSHong Zhang   Mat_MUMPS      *mumps;
1349ace3abfcSBarry Smith   PetscBool      isSeqAIJ;
13502877fffaSHong Zhang 
13512877fffaSHong Zhang   PetscFunctionBegin;
13522877fffaSHong Zhang   /* Create the factorization matrix */
1353251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr);
13542877fffaSHong Zhang   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
13552877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
13562877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1357bccb9932SShri Abhyankar   if (isSeqAIJ) {
13582877fffaSHong Zhang     ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
1359bccb9932SShri Abhyankar   } else {
1360bccb9932SShri Abhyankar     ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
1361bccb9932SShri Abhyankar   }
13622877fffaSHong Zhang 
136316ebf90aSShri Abhyankar   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
13642877fffaSHong Zhang   B->ops->view             = MatView_MUMPS;
136535bd34faSBarry Smith   B->ops->getinfo          = MatGetInfo_MUMPS;
136635bd34faSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
13675ccb76cbSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl_MUMPS",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
1368450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
1369450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
1370d5f3da31SBarry Smith     B->factortype = MAT_FACTOR_LU;
1371bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
1372bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
1373746480a1SHong Zhang     mumps->sym = 0;
1374dcd589f8SShri Abhyankar   } else {
137567877ebaSShri Abhyankar     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
1376450b117fSShri Abhyankar     B->factortype = MAT_FACTOR_CHOLESKY;
1377bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
1378bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
13796fdc2a6dSBarry Smith     if (A->spd_set && A->spd) mumps->sym = 1;
13806fdc2a6dSBarry Smith     else                      mumps->sym = 2;
1381450b117fSShri Abhyankar   }
13822877fffaSHong Zhang 
13832877fffaSHong Zhang   mumps->isAIJ        = PETSC_TRUE;
1384bf0cc555SLisandro Dalcin   mumps->Destroy      = B->ops->destroy;
13852877fffaSHong Zhang   B->ops->destroy     = MatDestroy_MUMPS;
13862877fffaSHong Zhang   B->spptr            = (void*)mumps;
1387f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1388746480a1SHong Zhang 
13892877fffaSHong Zhang   *F = B;
13902877fffaSHong Zhang   PetscFunctionReturn(0);
13912877fffaSHong Zhang }
13922877fffaSHong Zhang EXTERN_C_END
13932877fffaSHong Zhang 
1394bccb9932SShri Abhyankar 
13952877fffaSHong Zhang EXTERN_C_BEGIN
1396bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI SBAIJ matrices */
13972877fffaSHong Zhang #undef __FUNCT__
1398bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_sbaij_mumps"
1399bccb9932SShri Abhyankar PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
14002877fffaSHong Zhang {
14012877fffaSHong Zhang   Mat            B;
14022877fffaSHong Zhang   PetscErrorCode ierr;
14032877fffaSHong Zhang   Mat_MUMPS      *mumps;
1404ace3abfcSBarry Smith   PetscBool      isSeqSBAIJ;
14052877fffaSHong Zhang 
14062877fffaSHong Zhang   PetscFunctionBegin;
1407e7e72b3dSBarry Smith   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
1408bccb9932SShri Abhyankar   if (A->rmap->bs > 1) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with block size > 1 with MUMPS Cholesky, use AIJ matrix instead");
1409251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
14102877fffaSHong Zhang   /* Create the factorization matrix */
14112877fffaSHong Zhang   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
14122877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
14132877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
141416ebf90aSShri Abhyankar   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
1415bccb9932SShri Abhyankar   if (isSeqSBAIJ) {
1416bccb9932SShri Abhyankar     ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
141716ebf90aSShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
1418dcd589f8SShri Abhyankar   } else {
1419bccb9932SShri Abhyankar     ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
1420bccb9932SShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
1421bccb9932SShri Abhyankar   }
1422bccb9932SShri Abhyankar 
142367877ebaSShri Abhyankar   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
1424bccb9932SShri Abhyankar   B->ops->view                   = MatView_MUMPS;
1425bccb9932SShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
1426f250808bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl",MatMumpsSetIcntl);CHKERRQ(ierr);
1427f4762488SHong Zhang   B->factortype                  = MAT_FACTOR_CHOLESKY;
14286fdc2a6dSBarry Smith   if (A->spd_set && A->spd) mumps->sym = 1;
14296fdc2a6dSBarry Smith   else                      mumps->sym = 2;
1430a214ac2aSShri Abhyankar 
1431bccb9932SShri Abhyankar   mumps->isAIJ        = PETSC_FALSE;
1432bf0cc555SLisandro Dalcin   mumps->Destroy      = B->ops->destroy;
1433f3c0ef26SHong Zhang   B->ops->destroy     = MatDestroy_MUMPS;
14342877fffaSHong Zhang   B->spptr            = (void*)mumps;
1435f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1436746480a1SHong Zhang 
14372877fffaSHong Zhang   *F = B;
14382877fffaSHong Zhang   PetscFunctionReturn(0);
14392877fffaSHong Zhang }
14402877fffaSHong Zhang EXTERN_C_END
144197969023SHong Zhang 
1442450b117fSShri Abhyankar EXTERN_C_BEGIN
1443450b117fSShri Abhyankar #undef __FUNCT__
1444bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_baij_mumps"
1445bccb9932SShri Abhyankar PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
144667877ebaSShri Abhyankar {
144767877ebaSShri Abhyankar   Mat            B;
144867877ebaSShri Abhyankar   PetscErrorCode ierr;
144967877ebaSShri Abhyankar   Mat_MUMPS      *mumps;
1450ace3abfcSBarry Smith   PetscBool      isSeqBAIJ;
145167877ebaSShri Abhyankar 
145267877ebaSShri Abhyankar   PetscFunctionBegin;
145367877ebaSShri Abhyankar   /* Create the factorization matrix */
1454251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr);
145567877ebaSShri Abhyankar   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
145667877ebaSShri Abhyankar   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
145767877ebaSShri Abhyankar   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1458bccb9932SShri Abhyankar   if (isSeqBAIJ) {
145967877ebaSShri Abhyankar     ierr = MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,PETSC_NULL);CHKERRQ(ierr);
1460bccb9932SShri Abhyankar   } else {
146167877ebaSShri Abhyankar     ierr = MatMPIBAIJSetPreallocation(B,A->rmap->bs,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
1462bccb9932SShri Abhyankar   }
1463450b117fSShri Abhyankar 
146467877ebaSShri Abhyankar   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
1465450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
1466450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
1467450b117fSShri Abhyankar     B->factortype = MAT_FACTOR_LU;
1468bccb9932SShri Abhyankar     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
1469bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
1470746480a1SHong Zhang     mumps->sym = 0;
1471746480a1SHong Zhang   } else {
1472746480a1SHong Zhang     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n");
1473450b117fSShri Abhyankar   }
1474bccb9932SShri Abhyankar 
1475450b117fSShri Abhyankar   B->ops->view             = MatView_MUMPS;
1476450b117fSShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
14775ccb76cbSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl_MUMPS",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
1478450b117fSShri Abhyankar 
1479450b117fSShri Abhyankar   mumps->isAIJ        = PETSC_TRUE;
1480bf0cc555SLisandro Dalcin   mumps->Destroy      = B->ops->destroy;
1481450b117fSShri Abhyankar   B->ops->destroy     = MatDestroy_MUMPS;
1482450b117fSShri Abhyankar   B->spptr            = (void*)mumps;
1483f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1484746480a1SHong Zhang 
1485450b117fSShri Abhyankar   *F = B;
1486450b117fSShri Abhyankar   PetscFunctionReturn(0);
1487450b117fSShri Abhyankar }
1488450b117fSShri Abhyankar EXTERN_C_END
1489a214ac2aSShri Abhyankar 
1490