xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision 51d5961a1fd742257cce02e4a32c1ca9b3db276c)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
21c2a3de1SBarry Smith 
3397b6df1SKris Buschelman /*
4c2b5dc30SHong Zhang     Provides an interface to the MUMPS sparse solver
5397b6df1SKris Buschelman */
6*51d5961aSHong Zhang 
7*51d5961aSHong Zhang #include "../src/mat/impls/aij/mpi/mpiaij.h" /*I  "petscmat.h"  I*/
87c4f633dSBarry Smith #include "../src/mat/impls/sbaij/mpi/mpisbaij.h"
9397b6df1SKris Buschelman 
10397b6df1SKris Buschelman EXTERN_C_BEGIN
11397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
12397b6df1SKris Buschelman #include "zmumps_c.h"
13397b6df1SKris Buschelman #else
14397b6df1SKris Buschelman #include "dmumps_c.h"
15397b6df1SKris Buschelman #endif
16397b6df1SKris Buschelman EXTERN_C_END
17397b6df1SKris Buschelman #define JOB_INIT -1
183d472b54SHong Zhang #define JOB_FACTSYMBOLIC 1
193d472b54SHong Zhang #define JOB_FACTNUMERIC 2
203d472b54SHong Zhang #define JOB_SOLVE 3
21397b6df1SKris Buschelman #define JOB_END -2
223d472b54SHong Zhang 
233d472b54SHong Zhang 
24397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */
25397b6df1SKris Buschelman #define ICNTL(I) icntl[(I)-1]
26397b6df1SKris Buschelman #define CNTL(I) cntl[(I)-1]
27397b6df1SKris Buschelman #define INFOG(I) infog[(I)-1]
28a7aca84bSHong Zhang #define INFO(I) info[(I)-1]
29397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1]
30adc1d99fSHong Zhang #define RINFO(I) rinfo[(I)-1]
31397b6df1SKris Buschelman 
32397b6df1SKris Buschelman typedef struct {
33397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
34397b6df1SKris Buschelman   ZMUMPS_STRUC_C id;
35397b6df1SKris Buschelman #else
36397b6df1SKris Buschelman   DMUMPS_STRUC_C id;
37397b6df1SKris Buschelman #endif
38397b6df1SKris Buschelman   MatStructure   matstruc;
39c1490034SHong Zhang   PetscMPIInt    myid,size;
4016ebf90aSShri Abhyankar   PetscInt       *irn,*jcn,nz,sym,nSolve;
41397b6df1SKris Buschelman   PetscScalar    *val;
42397b6df1SKris Buschelman   MPI_Comm       comm_mumps;
43329ec9b3SHong Zhang   VecScatter     scat_rhs, scat_sol;
4464e6c443SBarry Smith   PetscBool      isAIJ,CleanUpMUMPS;
45329ec9b3SHong Zhang   Vec            b_seq,x_seq;
4667334b25SHong Zhang   PetscErrorCode (*MatDestroy)(Mat);
47bccb9932SShri Abhyankar   PetscErrorCode (*ConvertToTriples)(Mat, int, MatReuse, int*, int**, int**, PetscScalar**);
48f0c56d0fSKris Buschelman } Mat_MUMPS;
49f0c56d0fSKris Buschelman 
5009573ac7SBarry Smith extern PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);
51b24902e0SBarry Smith 
5267877ebaSShri Abhyankar 
5367877ebaSShri Abhyankar /* MatConvertToTriples_A_B */
5467877ebaSShri Abhyankar /*convert Petsc matrix to triples: row[nz], col[nz], val[nz] */
55397b6df1SKris Buschelman /*
56397b6df1SKris Buschelman   input:
5767877ebaSShri Abhyankar     A       - matrix in aij,baij or sbaij (bs=1) format
58397b6df1SKris Buschelman     shift   - 0: C style output triple; 1: Fortran style output triple.
59bccb9932SShri Abhyankar     reuse   - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple
60bccb9932SShri Abhyankar               MAT_REUSE_MATRIX:   only the values in v array are updated
61397b6df1SKris Buschelman   output:
62397b6df1SKris Buschelman     nnz     - dim of r, c, and v (number of local nonzero entries of A)
63397b6df1SKris Buschelman     r, c, v - row and col index, matrix values (matrix triples)
64397b6df1SKris Buschelman  */
6516ebf90aSShri Abhyankar 
6616ebf90aSShri Abhyankar #undef __FUNCT__
6716ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqaij"
68bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
69b24902e0SBarry Smith {
70185f6596SHong Zhang   const PetscInt   *ai,*aj,*ajj,M=A->rmap->n;
7167877ebaSShri Abhyankar   PetscInt         nz,rnz,i,j;
72dfbe8321SBarry Smith   PetscErrorCode   ierr;
73c1490034SHong Zhang   PetscInt         *row,*col;
7416ebf90aSShri Abhyankar   Mat_SeqAIJ       *aa=(Mat_SeqAIJ*)A->data;
75397b6df1SKris Buschelman 
76397b6df1SKris Buschelman   PetscFunctionBegin;
7716ebf90aSShri Abhyankar   *v=aa->a;
78bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX){
7916ebf90aSShri Abhyankar     nz = aa->nz; ai = aa->i; aj = aa->j;
8016ebf90aSShri Abhyankar     *nnz = nz;
81185f6596SHong Zhang     ierr = PetscMalloc(2*nz*sizeof(PetscInt), &row);CHKERRQ(ierr);
82185f6596SHong Zhang     col  = row + nz;
83185f6596SHong Zhang 
8416ebf90aSShri Abhyankar     nz = 0;
8516ebf90aSShri Abhyankar     for(i=0; i<M; i++) {
8616ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
8767877ebaSShri Abhyankar       ajj = aj + ai[i];
8867877ebaSShri Abhyankar       for(j=0; j<rnz; j++) {
8967877ebaSShri Abhyankar 	row[nz] = i+shift; col[nz++] = ajj[j] + shift;
9016ebf90aSShri Abhyankar       }
9116ebf90aSShri Abhyankar     }
9216ebf90aSShri Abhyankar     *r = row; *c = col;
9316ebf90aSShri Abhyankar   }
9416ebf90aSShri Abhyankar   PetscFunctionReturn(0);
9516ebf90aSShri Abhyankar }
96397b6df1SKris Buschelman 
9716ebf90aSShri Abhyankar #undef __FUNCT__
9867877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqbaij_seqaij"
99bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
10067877ebaSShri Abhyankar {
10167877ebaSShri Abhyankar   Mat_SeqBAIJ        *aa=(Mat_SeqBAIJ*)A->data;
10267877ebaSShri Abhyankar   const PetscInt     *ai,*aj,*ajj,bs=A->rmap->bs,bs2=aa->bs2,M=A->rmap->N/bs;
10367877ebaSShri Abhyankar   PetscInt           nz,idx=0,rnz,i,j,k,m,ii;
10467877ebaSShri Abhyankar   PetscErrorCode     ierr;
10567877ebaSShri Abhyankar   PetscInt           *row,*col;
10667877ebaSShri Abhyankar 
10767877ebaSShri Abhyankar   PetscFunctionBegin;
108cf3759fdSShri Abhyankar   *v = aa->a;
109bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX){
110cf3759fdSShri Abhyankar     ai = aa->i; aj = aa->j;
11167877ebaSShri Abhyankar     nz = bs2*aa->nz;
11267877ebaSShri Abhyankar     *nnz = nz;
113185f6596SHong Zhang     ierr = PetscMalloc(2*nz*sizeof(PetscInt), &row);CHKERRQ(ierr);
114185f6596SHong Zhang     col  = row + nz;
115185f6596SHong Zhang 
11667877ebaSShri Abhyankar     for(i=0; i<M; i++) {
11767877ebaSShri Abhyankar       ii = 0;
11867877ebaSShri Abhyankar       ajj = aj + ai[i];
11967877ebaSShri Abhyankar       rnz = ai[i+1] - ai[i];
12067877ebaSShri Abhyankar       for(k=0; k<rnz; k++) {
12167877ebaSShri Abhyankar 	for(j=0; j<bs; j++) {
12267877ebaSShri Abhyankar 	  for(m=0; m<bs; m++) {
12367877ebaSShri Abhyankar 	    row[idx]     = i*bs + m + shift;
124cf3759fdSShri Abhyankar 	    col[idx++]   = bs*(ajj[k]) + j + shift;
12567877ebaSShri Abhyankar 	  }
12667877ebaSShri Abhyankar 	}
12767877ebaSShri Abhyankar       }
12867877ebaSShri Abhyankar     }
129cf3759fdSShri Abhyankar     *r = row; *c = col;
13067877ebaSShri Abhyankar   }
13167877ebaSShri Abhyankar   PetscFunctionReturn(0);
13267877ebaSShri Abhyankar }
13367877ebaSShri Abhyankar 
13467877ebaSShri Abhyankar #undef __FUNCT__
13516ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqsbaij_seqsbaij"
136bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
13716ebf90aSShri Abhyankar {
13867877ebaSShri Abhyankar   const PetscInt   *ai, *aj,*ajj,M=A->rmap->n;
13967877ebaSShri Abhyankar   PetscInt         nz,rnz,i,j;
14016ebf90aSShri Abhyankar   PetscErrorCode   ierr;
14116ebf90aSShri Abhyankar   PetscInt         *row,*col;
14216ebf90aSShri Abhyankar   Mat_SeqSBAIJ     *aa=(Mat_SeqSBAIJ*)A->data;
14316ebf90aSShri Abhyankar 
14416ebf90aSShri Abhyankar   PetscFunctionBegin;
145bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX){
14616ebf90aSShri Abhyankar     nz = aa->nz;ai=aa->i; aj=aa->j;*v=aa->a;
14716ebf90aSShri Abhyankar     *nnz = nz;
148185f6596SHong Zhang     ierr = PetscMalloc(2*nz*sizeof(PetscInt), &row);CHKERRQ(ierr);
149185f6596SHong Zhang     col  = row + nz;
150185f6596SHong Zhang 
15116ebf90aSShri Abhyankar     nz = 0;
15216ebf90aSShri Abhyankar     for(i=0; i<M; i++) {
15316ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
15467877ebaSShri Abhyankar       ajj = aj + ai[i];
15567877ebaSShri Abhyankar       for(j=0; j<rnz; j++) {
15667877ebaSShri Abhyankar 	row[nz] = i+shift; col[nz++] = ajj[j] + shift;
15716ebf90aSShri Abhyankar       }
15816ebf90aSShri Abhyankar     }
15916ebf90aSShri Abhyankar     *r = row; *c = col;
16016ebf90aSShri Abhyankar   }
16116ebf90aSShri Abhyankar   PetscFunctionReturn(0);
16216ebf90aSShri Abhyankar }
16316ebf90aSShri Abhyankar 
16416ebf90aSShri Abhyankar #undef __FUNCT__
16516ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqsbaij"
166bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
16716ebf90aSShri Abhyankar {
16867877ebaSShri Abhyankar   const PetscInt     *ai,*aj,*ajj,*adiag,M=A->rmap->n;
16967877ebaSShri Abhyankar   PetscInt           nz,rnz,i,j;
17067877ebaSShri Abhyankar   const PetscScalar  *av,*v1;
17116ebf90aSShri Abhyankar   PetscScalar        *val;
17216ebf90aSShri Abhyankar   PetscErrorCode     ierr;
17316ebf90aSShri Abhyankar   PetscInt           *row,*col;
17416ebf90aSShri Abhyankar   Mat_SeqSBAIJ       *aa=(Mat_SeqSBAIJ*)A->data;
17516ebf90aSShri Abhyankar 
17616ebf90aSShri Abhyankar   PetscFunctionBegin;
17716ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j;av=aa->a;
17816ebf90aSShri Abhyankar   adiag=aa->diag;
179bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX){
18016ebf90aSShri Abhyankar     nz = M + (aa->nz-M)/2;
18116ebf90aSShri Abhyankar     *nnz = nz;
182185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
183185f6596SHong Zhang     col  = row + nz;
184185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
185185f6596SHong Zhang 
18616ebf90aSShri Abhyankar     nz = 0;
18716ebf90aSShri Abhyankar     for(i=0; i<M; i++) {
18816ebf90aSShri Abhyankar       rnz = ai[i+1] - adiag[i];
18967877ebaSShri Abhyankar       ajj  = aj + adiag[i];
190cf3759fdSShri Abhyankar       v1   = av + adiag[i];
19167877ebaSShri Abhyankar       for(j=0; j<rnz; j++) {
19267877ebaSShri Abhyankar 	row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
19316ebf90aSShri Abhyankar       }
19416ebf90aSShri Abhyankar     }
19516ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
196397b6df1SKris Buschelman   } else {
19716ebf90aSShri Abhyankar     nz = 0; val = *v;
19816ebf90aSShri Abhyankar     for(i=0; i <M; i++) {
19916ebf90aSShri Abhyankar       rnz = ai[i+1] - adiag[i];
20067877ebaSShri Abhyankar       ajj = aj + adiag[i];
20167877ebaSShri Abhyankar       v1  = av + adiag[i];
20267877ebaSShri Abhyankar       for(j=0; j<rnz; j++) {
20367877ebaSShri Abhyankar 	val[nz++] = v1[j];
20416ebf90aSShri Abhyankar       }
20516ebf90aSShri Abhyankar     }
20616ebf90aSShri Abhyankar   }
20716ebf90aSShri Abhyankar   PetscFunctionReturn(0);
20816ebf90aSShri Abhyankar }
20916ebf90aSShri Abhyankar 
21016ebf90aSShri Abhyankar #undef __FUNCT__
21116ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpisbaij_mpisbaij"
212bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
21316ebf90aSShri Abhyankar {
21416ebf90aSShri Abhyankar   const PetscInt     *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
21516ebf90aSShri Abhyankar   PetscErrorCode     ierr;
21616ebf90aSShri Abhyankar   PetscInt           rstart,nz,i,j,jj,irow,countA,countB;
21716ebf90aSShri Abhyankar   PetscInt           *row,*col;
21816ebf90aSShri Abhyankar   const PetscScalar  *av, *bv,*v1,*v2;
21916ebf90aSShri Abhyankar   PetscScalar        *val;
220397b6df1SKris Buschelman   Mat_MPISBAIJ       *mat =  (Mat_MPISBAIJ*)A->data;
221397b6df1SKris Buschelman   Mat_SeqSBAIJ       *aa=(Mat_SeqSBAIJ*)(mat->A)->data;
222397b6df1SKris Buschelman   Mat_SeqBAIJ        *bb=(Mat_SeqBAIJ*)(mat->B)->data;
22316ebf90aSShri Abhyankar 
22416ebf90aSShri Abhyankar   PetscFunctionBegin;
225d0f46423SBarry Smith   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
226397b6df1SKris Buschelman   garray = mat->garray;
227397b6df1SKris Buschelman   av=aa->a; bv=bb->a;
228397b6df1SKris Buschelman 
229bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX){
23016ebf90aSShri Abhyankar     nz = aa->nz + bb->nz;
23116ebf90aSShri Abhyankar     *nnz = nz;
232185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
233185f6596SHong Zhang     col  = row + nz;
234185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
235185f6596SHong Zhang 
236397b6df1SKris Buschelman     *r = row; *c = col; *v = val;
237397b6df1SKris Buschelman   } else {
238397b6df1SKris Buschelman     row = *r; col = *c; val = *v;
239397b6df1SKris Buschelman   }
240397b6df1SKris Buschelman 
241028e57e8SHong Zhang   jj = 0; irow = rstart;
242397b6df1SKris Buschelman   for ( i=0; i<m; i++ ) {
243397b6df1SKris Buschelman     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
244397b6df1SKris Buschelman     countA = ai[i+1] - ai[i];
245397b6df1SKris Buschelman     countB = bi[i+1] - bi[i];
246397b6df1SKris Buschelman     bjj    = bj + bi[i];
24716ebf90aSShri Abhyankar     v1     = av + ai[i];
24816ebf90aSShri Abhyankar     v2     = bv + bi[i];
249397b6df1SKris Buschelman 
250397b6df1SKris Buschelman     /* A-part */
251397b6df1SKris Buschelman     for (j=0; j<countA; j++){
252bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
253397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
254397b6df1SKris Buschelman       }
25516ebf90aSShri Abhyankar       val[jj++] = v1[j];
256397b6df1SKris Buschelman     }
25716ebf90aSShri Abhyankar 
25816ebf90aSShri Abhyankar     /* B-part */
25916ebf90aSShri Abhyankar     for(j=0; j < countB; j++){
260bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
261397b6df1SKris Buschelman 	row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
262397b6df1SKris Buschelman       }
26316ebf90aSShri Abhyankar       val[jj++] = v2[j];
26416ebf90aSShri Abhyankar     }
26516ebf90aSShri Abhyankar     irow++;
26616ebf90aSShri Abhyankar   }
26716ebf90aSShri Abhyankar   PetscFunctionReturn(0);
26816ebf90aSShri Abhyankar }
26916ebf90aSShri Abhyankar 
27016ebf90aSShri Abhyankar #undef __FUNCT__
27116ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpiaij"
272bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
27316ebf90aSShri Abhyankar {
27416ebf90aSShri Abhyankar   const PetscInt     *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
27516ebf90aSShri Abhyankar   PetscErrorCode     ierr;
27616ebf90aSShri Abhyankar   PetscInt           rstart,nz,i,j,jj,irow,countA,countB;
27716ebf90aSShri Abhyankar   PetscInt           *row,*col;
27816ebf90aSShri Abhyankar   const PetscScalar  *av, *bv,*v1,*v2;
27916ebf90aSShri Abhyankar   PetscScalar        *val;
28016ebf90aSShri Abhyankar   Mat_MPIAIJ         *mat =  (Mat_MPIAIJ*)A->data;
28116ebf90aSShri Abhyankar   Mat_SeqAIJ         *aa=(Mat_SeqAIJ*)(mat->A)->data;
28216ebf90aSShri Abhyankar   Mat_SeqAIJ         *bb=(Mat_SeqAIJ*)(mat->B)->data;
28316ebf90aSShri Abhyankar 
28416ebf90aSShri Abhyankar   PetscFunctionBegin;
28516ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
28616ebf90aSShri Abhyankar   garray = mat->garray;
28716ebf90aSShri Abhyankar   av=aa->a; bv=bb->a;
28816ebf90aSShri Abhyankar 
289bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX){
29016ebf90aSShri Abhyankar     nz = aa->nz + bb->nz;
29116ebf90aSShri Abhyankar     *nnz = nz;
292185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
293185f6596SHong Zhang     col  = row + nz;
294185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
295185f6596SHong Zhang 
29616ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
29716ebf90aSShri Abhyankar   } else {
29816ebf90aSShri Abhyankar     row = *r; col = *c; val = *v;
29916ebf90aSShri Abhyankar   }
30016ebf90aSShri Abhyankar 
30116ebf90aSShri Abhyankar   jj = 0; irow = rstart;
30216ebf90aSShri Abhyankar   for ( i=0; i<m; i++ ) {
30316ebf90aSShri Abhyankar     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
30416ebf90aSShri Abhyankar     countA = ai[i+1] - ai[i];
30516ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
30616ebf90aSShri Abhyankar     bjj    = bj + bi[i];
30716ebf90aSShri Abhyankar     v1     = av + ai[i];
30816ebf90aSShri Abhyankar     v2     = bv + bi[i];
30916ebf90aSShri Abhyankar 
31016ebf90aSShri Abhyankar     /* A-part */
31116ebf90aSShri Abhyankar     for (j=0; j<countA; j++){
312bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX){
31316ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
31416ebf90aSShri Abhyankar       }
31516ebf90aSShri Abhyankar       val[jj++] = v1[j];
31616ebf90aSShri Abhyankar     }
31716ebf90aSShri Abhyankar 
31816ebf90aSShri Abhyankar     /* B-part */
31916ebf90aSShri Abhyankar     for(j=0; j < countB; j++){
320bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX){
32116ebf90aSShri Abhyankar 	row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
32216ebf90aSShri Abhyankar       }
32316ebf90aSShri Abhyankar       val[jj++] = v2[j];
32416ebf90aSShri Abhyankar     }
32516ebf90aSShri Abhyankar     irow++;
32616ebf90aSShri Abhyankar   }
32716ebf90aSShri Abhyankar   PetscFunctionReturn(0);
32816ebf90aSShri Abhyankar }
32916ebf90aSShri Abhyankar 
33016ebf90aSShri Abhyankar #undef __FUNCT__
33167877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpibaij_mpiaij"
332bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
33367877ebaSShri Abhyankar {
33467877ebaSShri Abhyankar   Mat_MPIBAIJ        *mat =  (Mat_MPIBAIJ*)A->data;
33567877ebaSShri Abhyankar   Mat_SeqBAIJ        *aa=(Mat_SeqBAIJ*)(mat->A)->data;
33667877ebaSShri Abhyankar   Mat_SeqBAIJ        *bb=(Mat_SeqBAIJ*)(mat->B)->data;
33767877ebaSShri Abhyankar   const PetscInt     *ai = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj;
33867877ebaSShri Abhyankar   const PetscInt     *garray = mat->garray,mbs=mat->mbs,rstartbs=mat->rstartbs;
33967877ebaSShri Abhyankar   const PetscInt     bs = A->rmap->bs,bs2=mat->bs2;
34067877ebaSShri Abhyankar   PetscErrorCode     ierr;
34167877ebaSShri Abhyankar   PetscInt           nz,i,j,k,n,jj,irow,countA,countB,idx;
34267877ebaSShri Abhyankar   PetscInt           *row,*col;
34367877ebaSShri Abhyankar   const PetscScalar  *av=aa->a, *bv=bb->a,*v1,*v2;
34467877ebaSShri Abhyankar   PetscScalar        *val;
34567877ebaSShri Abhyankar 
34667877ebaSShri Abhyankar   PetscFunctionBegin;
34767877ebaSShri Abhyankar 
348bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
34967877ebaSShri Abhyankar     nz = bs2*(aa->nz + bb->nz);
35067877ebaSShri Abhyankar     *nnz = nz;
351185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
352185f6596SHong Zhang     col  = row + nz;
353185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
354185f6596SHong Zhang 
35567877ebaSShri Abhyankar     *r = row; *c = col; *v = val;
35667877ebaSShri Abhyankar   } else {
35767877ebaSShri Abhyankar     row = *r; col = *c; val = *v;
35867877ebaSShri Abhyankar   }
35967877ebaSShri Abhyankar 
36067877ebaSShri Abhyankar   jj = 0; irow = rstartbs;
36167877ebaSShri Abhyankar   for ( i=0; i<mbs; i++ ) {
36267877ebaSShri Abhyankar     countA = ai[i+1] - ai[i];
36367877ebaSShri Abhyankar     countB = bi[i+1] - bi[i];
36467877ebaSShri Abhyankar     ajj    = aj + ai[i];
36567877ebaSShri Abhyankar     bjj    = bj + bi[i];
36667877ebaSShri Abhyankar     v1     = av + bs2*ai[i];
36767877ebaSShri Abhyankar     v2     = bv + bs2*bi[i];
36867877ebaSShri Abhyankar 
36967877ebaSShri Abhyankar     idx = 0;
37067877ebaSShri Abhyankar     /* A-part */
37167877ebaSShri Abhyankar     for (k=0; k<countA; k++){
37267877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
37367877ebaSShri Abhyankar 	for (n=0; n<bs; n++) {
374bccb9932SShri Abhyankar 	  if (reuse == MAT_INITIAL_MATRIX){
37567877ebaSShri Abhyankar 	    row[jj] = bs*irow + n + shift;
37667877ebaSShri Abhyankar 	    col[jj] = bs*(rstartbs + ajj[k]) + j + shift;
37767877ebaSShri Abhyankar 	  }
37867877ebaSShri Abhyankar 	  val[jj++] = v1[idx++];
37967877ebaSShri Abhyankar 	}
38067877ebaSShri Abhyankar       }
38167877ebaSShri Abhyankar     }
38267877ebaSShri Abhyankar 
38367877ebaSShri Abhyankar     idx = 0;
38467877ebaSShri Abhyankar     /* B-part */
38567877ebaSShri Abhyankar     for(k=0; k<countB; k++){
38667877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
38767877ebaSShri Abhyankar 	for (n=0; n<bs; n++) {
388bccb9932SShri Abhyankar 	  if (reuse == MAT_INITIAL_MATRIX){
38967877ebaSShri Abhyankar 	    row[jj] = bs*irow + n + shift;
39067877ebaSShri Abhyankar 	    col[jj] = bs*(garray[bjj[k]]) + j + shift;
39167877ebaSShri Abhyankar 	  }
39267877ebaSShri Abhyankar 	  val[jj++] = bv[idx++];
39367877ebaSShri Abhyankar 	}
39467877ebaSShri Abhyankar       }
39567877ebaSShri Abhyankar     }
39667877ebaSShri Abhyankar     irow++;
39767877ebaSShri Abhyankar   }
39867877ebaSShri Abhyankar   PetscFunctionReturn(0);
39967877ebaSShri Abhyankar }
40067877ebaSShri Abhyankar 
40167877ebaSShri Abhyankar #undef __FUNCT__
40216ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpisbaij"
403bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
40416ebf90aSShri Abhyankar {
40516ebf90aSShri Abhyankar   const PetscInt     *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
40616ebf90aSShri Abhyankar   PetscErrorCode     ierr;
40716ebf90aSShri Abhyankar   PetscInt           rstart,nz,nza,nzb_low,i,j,jj,irow,countA,countB;
40816ebf90aSShri Abhyankar   PetscInt           *row,*col;
40916ebf90aSShri Abhyankar   const PetscScalar  *av, *bv,*v1,*v2;
41016ebf90aSShri Abhyankar   PetscScalar        *val;
41116ebf90aSShri Abhyankar   Mat_MPIAIJ         *mat =  (Mat_MPIAIJ*)A->data;
41216ebf90aSShri Abhyankar   Mat_SeqAIJ         *aa=(Mat_SeqAIJ*)(mat->A)->data;
41316ebf90aSShri Abhyankar   Mat_SeqAIJ         *bb=(Mat_SeqAIJ*)(mat->B)->data;
41416ebf90aSShri Abhyankar 
41516ebf90aSShri Abhyankar   PetscFunctionBegin;
41616ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j; adiag=aa->diag;
41716ebf90aSShri Abhyankar   bi=bb->i; bj=bb->j; garray = mat->garray;
41816ebf90aSShri Abhyankar   av=aa->a; bv=bb->a;
41916ebf90aSShri Abhyankar   rstart = A->rmap->rstart;
42016ebf90aSShri Abhyankar 
421bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
42216ebf90aSShri Abhyankar     nza = 0;nzb_low = 0;
42316ebf90aSShri Abhyankar     for(i=0; i<m; i++){
42416ebf90aSShri Abhyankar       nza     = nza + (ai[i+1] - adiag[i]);
42516ebf90aSShri Abhyankar       countB  = bi[i+1] - bi[i];
42616ebf90aSShri Abhyankar       bjj     = bj + bi[i];
42716ebf90aSShri Abhyankar 
42816ebf90aSShri Abhyankar       j = 0;
42916ebf90aSShri Abhyankar       while(garray[bjj[j]] < rstart) {
43016ebf90aSShri Abhyankar 	if(j == countB) break;
43116ebf90aSShri Abhyankar 	j++;nzb_low++;
43216ebf90aSShri Abhyankar       }
43316ebf90aSShri Abhyankar     }
43416ebf90aSShri Abhyankar     /* Total nz = nz for the upper triangular A part + nz for the 2nd B part */
43516ebf90aSShri Abhyankar     nz = nza + (bb->nz - nzb_low);
43616ebf90aSShri Abhyankar     *nnz = nz;
437185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
438185f6596SHong Zhang     col  = row + nz;
439185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
440185f6596SHong Zhang 
44116ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
44216ebf90aSShri Abhyankar   } else {
44316ebf90aSShri Abhyankar     row = *r; col = *c; val = *v;
44416ebf90aSShri Abhyankar   }
44516ebf90aSShri Abhyankar 
44616ebf90aSShri Abhyankar   jj = 0; irow = rstart;
44716ebf90aSShri Abhyankar   for ( i=0; i<m; i++ ) {
44816ebf90aSShri Abhyankar     ajj    = aj + adiag[i];                 /* ptr to the beginning of the diagonal of this row */
44916ebf90aSShri Abhyankar     v1     = av + adiag[i];
45016ebf90aSShri Abhyankar     countA = ai[i+1] - adiag[i];
45116ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
45216ebf90aSShri Abhyankar     bjj    = bj + bi[i];
45316ebf90aSShri Abhyankar     v2     = bv + bi[i];
45416ebf90aSShri Abhyankar 
45516ebf90aSShri Abhyankar      /* A-part */
45616ebf90aSShri Abhyankar     for (j=0; j<countA; j++){
457bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
45816ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
45916ebf90aSShri Abhyankar       }
46016ebf90aSShri Abhyankar       val[jj++] = v1[j];
46116ebf90aSShri Abhyankar     }
46216ebf90aSShri Abhyankar 
46316ebf90aSShri Abhyankar     /* B-part */
46416ebf90aSShri Abhyankar     for(j=0; j < countB; j++){
46516ebf90aSShri Abhyankar       if (garray[bjj[j]] > rstart) {
466bccb9932SShri Abhyankar 	if (reuse == MAT_INITIAL_MATRIX) {
46716ebf90aSShri Abhyankar 	  row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
46816ebf90aSShri Abhyankar 	}
46916ebf90aSShri Abhyankar 	val[jj++] = v2[j];
47016ebf90aSShri Abhyankar       }
471397b6df1SKris Buschelman     }
472397b6df1SKris Buschelman     irow++;
473397b6df1SKris Buschelman   }
474397b6df1SKris Buschelman   PetscFunctionReturn(0);
475397b6df1SKris Buschelman }
476397b6df1SKris Buschelman 
477397b6df1SKris Buschelman #undef __FUNCT__
4783924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS"
479dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A)
480dfbe8321SBarry Smith {
481f0c56d0fSKris Buschelman   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
482dfbe8321SBarry Smith   PetscErrorCode ierr;
483b24902e0SBarry Smith 
484397b6df1SKris Buschelman   PetscFunctionBegin;
485397b6df1SKris Buschelman   if (lu->CleanUpMUMPS) {
486397b6df1SKris Buschelman     /* Terminate instance, deallocate memories */
4878fa425b9SSatish Balay     if (lu->id.sol_loc){ierr = PetscFree2(lu->id.sol_loc,lu->id.isol_loc);CHKERRQ(ierr);}
4888fa425b9SSatish Balay     if (lu->scat_rhs){ierr = VecScatterDestroy(lu->scat_rhs);CHKERRQ(ierr);}
4898fa425b9SSatish Balay     if (lu->b_seq) {ierr = VecDestroy(lu->b_seq);CHKERRQ(ierr);}
4902750af12SHong Zhang     if (lu->nSolve && lu->scat_sol){ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);}
4912750af12SHong Zhang     if (lu->nSolve && lu->x_seq){ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);}
4927c93a85dSSatish Balay 
493185f6596SHong Zhang     ierr = PetscFree(lu->irn);CHKERRQ(ierr);
494397b6df1SKris Buschelman     lu->id.job=JOB_END;
495397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
496397b6df1SKris Buschelman     zmumps_c(&lu->id);
497397b6df1SKris Buschelman #else
498397b6df1SKris Buschelman     dmumps_c(&lu->id);
499397b6df1SKris Buschelman #endif
500397b6df1SKris Buschelman     ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr);
501397b6df1SKris Buschelman   }
50297969023SHong Zhang   /* clear composed functions */
50397969023SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatFactorGetSolverPackage_C","",PETSC_NULL);CHKERRQ(ierr);
504f250808bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMumpsSetIcntl_C","",PETSC_NULL);CHKERRQ(ierr);
50567334b25SHong Zhang   ierr = (lu->MatDestroy)(A);CHKERRQ(ierr);
506397b6df1SKris Buschelman   PetscFunctionReturn(0);
507397b6df1SKris Buschelman }
508397b6df1SKris Buschelman 
509397b6df1SKris Buschelman #undef __FUNCT__
510f6c57405SHong Zhang #define __FUNCT__ "MatSolve_MUMPS"
511b24902e0SBarry Smith PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x)
512b24902e0SBarry Smith {
513f0c56d0fSKris Buschelman   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
514d54de34fSKris Buschelman   PetscScalar    *array;
51567877ebaSShri Abhyankar   Vec            b_seq;
516329ec9b3SHong Zhang   IS             is_iden,is_petsc;
517dfbe8321SBarry Smith   PetscErrorCode ierr;
518329ec9b3SHong Zhang   PetscInt       i;
519397b6df1SKris Buschelman 
520397b6df1SKris Buschelman   PetscFunctionBegin;
521329ec9b3SHong Zhang   lu->id.nrhs = 1;
52267877ebaSShri Abhyankar   b_seq = lu->b_seq;
523397b6df1SKris Buschelman   if (lu->size > 1){
524329ec9b3SHong Zhang     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
52567877ebaSShri Abhyankar     ierr = VecScatterBegin(lu->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
52667877ebaSShri Abhyankar     ierr = VecScatterEnd(lu->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
52767877ebaSShri Abhyankar     if (!lu->myid) {ierr = VecGetArray(b_seq,&array);CHKERRQ(ierr);}
528397b6df1SKris Buschelman   } else {  /* size == 1 */
529397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
530397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
531397b6df1SKris Buschelman   }
532397b6df1SKris Buschelman   if (!lu->myid) { /* define rhs on the host */
5338278f211SHong Zhang     lu->id.nrhs = 1;
534397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
535397b6df1SKris Buschelman     lu->id.rhs = (mumps_double_complex*)array;
536397b6df1SKris Buschelman #else
537397b6df1SKris Buschelman     lu->id.rhs = array;
538397b6df1SKris Buschelman #endif
539397b6df1SKris Buschelman   }
540397b6df1SKris Buschelman 
541397b6df1SKris Buschelman   /* solve phase */
542329ec9b3SHong Zhang   /*-------------*/
5433d472b54SHong Zhang   lu->id.job = JOB_SOLVE;
544397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
545397b6df1SKris Buschelman   zmumps_c(&lu->id);
546397b6df1SKris Buschelman #else
547397b6df1SKris Buschelman   dmumps_c(&lu->id);
548397b6df1SKris Buschelman #endif
54965e19b50SBarry 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));
550397b6df1SKris Buschelman 
551329ec9b3SHong Zhang   if (lu->size > 1) { /* convert mumps distributed solution to petsc mpi x */
552329ec9b3SHong Zhang     if (!lu->nSolve){ /* create scatter scat_sol */
553329ec9b3SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,lu->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */
554329ec9b3SHong Zhang       for (i=0; i<lu->id.lsol_loc; i++){
555329ec9b3SHong Zhang         lu->id.isol_loc[i] -= 1; /* change Fortran style to C style */
556397b6df1SKris Buschelman       }
55770b3c8c7SBarry Smith       ierr = ISCreateGeneral(PETSC_COMM_SELF,lu->id.lsol_loc,lu->id.isol_loc,PETSC_COPY_VALUES,&is_petsc);CHKERRQ(ierr);  /* to */
558329ec9b3SHong Zhang       ierr = VecScatterCreate(lu->x_seq,is_iden,x,is_petsc,&lu->scat_sol);CHKERRQ(ierr);
559329ec9b3SHong Zhang       ierr = ISDestroy(is_iden);CHKERRQ(ierr);
560329ec9b3SHong Zhang       ierr = ISDestroy(is_petsc);CHKERRQ(ierr);
561397b6df1SKris Buschelman     }
562ca9f406cSSatish Balay     ierr = VecScatterBegin(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
563ca9f406cSSatish Balay     ierr = VecScatterEnd(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
564329ec9b3SHong Zhang   }
565329ec9b3SHong Zhang   lu->nSolve++;
566397b6df1SKris Buschelman   PetscFunctionReturn(0);
567397b6df1SKris Buschelman }
568397b6df1SKris Buschelman 
569*51d5961aSHong Zhang #undef __FUNCT__
570*51d5961aSHong Zhang #define __FUNCT__ "MatSolveTranspose_MUMPS"
571*51d5961aSHong Zhang PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x)
572*51d5961aSHong Zhang {
573*51d5961aSHong Zhang   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
574*51d5961aSHong Zhang   PetscErrorCode ierr;
575*51d5961aSHong Zhang 
576*51d5961aSHong Zhang   PetscFunctionBegin;
577*51d5961aSHong Zhang   lu->id.ICNTL(9) = 0;
578*51d5961aSHong Zhang   ierr = MatSolve_MUMPS(A,b,x);
579*51d5961aSHong Zhang   lu->id.ICNTL(9) = 1;
580*51d5961aSHong Zhang   PetscFunctionReturn(0);
581*51d5961aSHong Zhang }
582*51d5961aSHong Zhang 
583ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX)
584a58c3f20SHong Zhang /*
585a58c3f20SHong Zhang   input:
586a58c3f20SHong Zhang    F:        numeric factor
587a58c3f20SHong Zhang   output:
588a58c3f20SHong Zhang    nneg:     total number of negative pivots
589a58c3f20SHong Zhang    nzero:    0
590a58c3f20SHong Zhang    npos:     (global dimension of F) - nneg
591a58c3f20SHong Zhang */
592a58c3f20SHong Zhang 
593a58c3f20SHong Zhang #undef __FUNCT__
594a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
595dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
596a58c3f20SHong Zhang {
597a58c3f20SHong Zhang   Mat_MUMPS      *lu =(Mat_MUMPS*)F->spptr;
598dfbe8321SBarry Smith   PetscErrorCode ierr;
599c1490034SHong Zhang   PetscMPIInt    size;
600a58c3f20SHong Zhang 
601a58c3f20SHong Zhang   PetscFunctionBegin;
6027adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)F)->comm,&size);CHKERRQ(ierr);
603bcb30aebSHong 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 */
60465e19b50SBarry 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));
605a58c3f20SHong Zhang   if (nneg){
606a58c3f20SHong Zhang     if (!lu->myid){
607a58c3f20SHong Zhang       *nneg = lu->id.INFOG(12);
608a58c3f20SHong Zhang     }
609bcb30aebSHong Zhang     ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr);
610a58c3f20SHong Zhang   }
611a58c3f20SHong Zhang   if (nzero) *nzero = 0;
612d0f46423SBarry Smith   if (npos)  *npos  = F->rmap->N - (*nneg);
613a58c3f20SHong Zhang   PetscFunctionReturn(0);
614a58c3f20SHong Zhang }
615ace3df97SHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */
616a58c3f20SHong Zhang 
617397b6df1SKris Buschelman #undef __FUNCT__
618f6c57405SHong Zhang #define __FUNCT__ "MatFactorNumeric_MUMPS"
6190481f469SBarry Smith PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
620af281ebdSHong Zhang {
621dcd589f8SShri Abhyankar   Mat_MUMPS       *lu =(Mat_MUMPS*)(F)->spptr;
6226849ba73SBarry Smith   PetscErrorCode  ierr;
623bccb9932SShri Abhyankar   MatReuse        reuse;
624e09efc27SHong Zhang   Mat             F_diag;
625ace3abfcSBarry Smith   PetscBool       isMPIAIJ;
626397b6df1SKris Buschelman 
627397b6df1SKris Buschelman   PetscFunctionBegin;
628bccb9932SShri Abhyankar   reuse = MAT_REUSE_MATRIX;
629bccb9932SShri Abhyankar   ierr = (*lu->ConvertToTriples)(A, 1, reuse, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
630397b6df1SKris Buschelman 
631397b6df1SKris Buschelman   /* numerical factorization phase */
632329ec9b3SHong Zhang   /*-------------------------------*/
6333d472b54SHong Zhang   lu->id.job = JOB_FACTNUMERIC;
634958c9bccSBarry Smith   if(!lu->id.ICNTL(18)) {
635a7aca84bSHong Zhang     if (!lu->myid) {
636397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
637397b6df1SKris Buschelman       lu->id.a = (mumps_double_complex*)lu->val;
638397b6df1SKris Buschelman #else
639397b6df1SKris Buschelman       lu->id.a = lu->val;
640397b6df1SKris Buschelman #endif
641397b6df1SKris Buschelman     }
642397b6df1SKris Buschelman   } else {
643397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
644397b6df1SKris Buschelman     lu->id.a_loc = (mumps_double_complex*)lu->val;
645397b6df1SKris Buschelman #else
646397b6df1SKris Buschelman     lu->id.a_loc = lu->val;
647397b6df1SKris Buschelman #endif
648397b6df1SKris Buschelman   }
649397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
650397b6df1SKris Buschelman   zmumps_c(&lu->id);
651397b6df1SKris Buschelman #else
652397b6df1SKris Buschelman   dmumps_c(&lu->id);
653397b6df1SKris Buschelman #endif
654397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
65565e19b50SBarry 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));
65665e19b50SBarry 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));
657397b6df1SKris Buschelman   }
65865e19b50SBarry 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));
659397b6df1SKris Buschelman 
6608ada1bb4SHong Zhang   if (lu->size > 1){
66116ebf90aSShri Abhyankar     ierr = PetscTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr);
662a214ac2aSShri Abhyankar     if(isMPIAIJ) {
663dcd589f8SShri Abhyankar       F_diag = ((Mat_MPIAIJ *)(F)->data)->A;
664e09efc27SHong Zhang     } else {
665dcd589f8SShri Abhyankar       F_diag = ((Mat_MPISBAIJ *)(F)->data)->A;
666e09efc27SHong Zhang     }
667e09efc27SHong Zhang     F_diag->assembled = PETSC_TRUE;
668329ec9b3SHong Zhang     if (lu->nSolve){
669329ec9b3SHong Zhang       ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);
6700e83c824SBarry Smith       ierr = PetscFree2(lu->id.sol_loc,lu->id.isol_loc);CHKERRQ(ierr);
671329ec9b3SHong Zhang       ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);
672329ec9b3SHong Zhang     }
6738ada1bb4SHong Zhang   }
674dcd589f8SShri Abhyankar   (F)->assembled   = PETSC_TRUE;
675397b6df1SKris Buschelman   lu->matstruc     = SAME_NONZERO_PATTERN;
676ace87b0dSHong Zhang   lu->CleanUpMUMPS = PETSC_TRUE;
677329ec9b3SHong Zhang   lu->nSolve       = 0;
67867877ebaSShri Abhyankar 
67967877ebaSShri Abhyankar   if (lu->size > 1){
68067877ebaSShri Abhyankar     /* distributed solution */
68167877ebaSShri Abhyankar     lu->id.ICNTL(21) = 1;
68267877ebaSShri Abhyankar     if (!lu->nSolve){
68367877ebaSShri Abhyankar       /* Create x_seq=sol_loc for repeated use */
68467877ebaSShri Abhyankar       PetscInt    lsol_loc;
68567877ebaSShri Abhyankar       PetscScalar *sol_loc;
68667877ebaSShri Abhyankar       lsol_loc = lu->id.INFO(23); /* length of sol_loc */
68767877ebaSShri Abhyankar       ierr = PetscMalloc2(lsol_loc,PetscScalar,&sol_loc,lsol_loc,PetscInt,&lu->id.isol_loc);CHKERRQ(ierr);
68867877ebaSShri Abhyankar       lu->id.lsol_loc = lsol_loc;
68967877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
69067877ebaSShri Abhyankar       lu->id.sol_loc  = (mumps_double_complex*)sol_loc;
69167877ebaSShri Abhyankar #else
69267877ebaSShri Abhyankar       lu->id.sol_loc  = sol_loc;
69367877ebaSShri Abhyankar #endif
69467877ebaSShri Abhyankar       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,lsol_loc,sol_loc,&lu->x_seq);CHKERRQ(ierr);
69567877ebaSShri Abhyankar     }
69667877ebaSShri Abhyankar   }
697397b6df1SKris Buschelman   PetscFunctionReturn(0);
698397b6df1SKris Buschelman }
699397b6df1SKris Buschelman 
700dcd589f8SShri Abhyankar #undef __FUNCT__
701dcd589f8SShri Abhyankar #define __FUNCT__ "PetscSetMUMPSOptions"
702dcd589f8SShri Abhyankar PetscErrorCode PetscSetMUMPSOptions(Mat F, Mat A)
703dcd589f8SShri Abhyankar {
704dcd589f8SShri Abhyankar   Mat_MUMPS        *lu = (Mat_MUMPS*)F->spptr;
705dcd589f8SShri Abhyankar   PetscErrorCode   ierr;
706dcd589f8SShri Abhyankar   PetscInt         icntl;
707ace3abfcSBarry Smith   PetscBool        flg;
708dcd589f8SShri Abhyankar 
709dcd589f8SShri Abhyankar   PetscFunctionBegin;
710dcd589f8SShri Abhyankar   ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
711dcd589f8SShri Abhyankar   if (lu->size == 1){
712dcd589f8SShri Abhyankar     lu->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
713dcd589f8SShri Abhyankar   } else {
714dcd589f8SShri Abhyankar     lu->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
715dcd589f8SShri Abhyankar   }
716dcd589f8SShri Abhyankar 
717dcd589f8SShri Abhyankar   icntl=-1;
718dcd589f8SShri Abhyankar   lu->id.ICNTL(4) = 0;  /* level of printing; overwrite mumps default ICNTL(4)=2 */
719dcd589f8SShri Abhyankar   ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
720dcd589f8SShri Abhyankar   if ((flg && icntl > 0) || PetscLogPrintInfo) {
721dcd589f8SShri Abhyankar     lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */
722dcd589f8SShri Abhyankar   } else { /* no output */
723dcd589f8SShri Abhyankar     lu->id.ICNTL(1) = 0;  /* error message, default= 6 */
724dcd589f8SShri Abhyankar     lu->id.ICNTL(2) = 0;  /* output stream for diagnostic printing, statistics, and warning. default=0 */
725dcd589f8SShri Abhyankar     lu->id.ICNTL(3) = 0; /* output stream for global information, default=6 */
726dcd589f8SShri Abhyankar   }
727dcd589f8SShri Abhyankar   ierr = PetscOptionsInt("-mat_mumps_icntl_6","ICNTL(6): column permutation and/or scaling to get a zero-free diagonal (0 to 7)","None",lu->id.ICNTL(6),&lu->id.ICNTL(6),PETSC_NULL);CHKERRQ(ierr);
728dcd589f8SShri Abhyankar   icntl=-1;
729292fb18eSBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): sequential matrix ordering (0 to 7) 3 = Scotch, 5 = Metis","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
730dcd589f8SShri Abhyankar   if (flg) {
731dcd589f8SShri Abhyankar     if (icntl== 1){
732e32f2f54SBarry 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");
733dcd589f8SShri Abhyankar     } else {
734dcd589f8SShri Abhyankar       lu->id.ICNTL(7) = icntl;
735dcd589f8SShri Abhyankar     }
736dcd589f8SShri Abhyankar   }
737dcd589f8SShri Abhyankar   ierr = PetscOptionsInt("-mat_mumps_icntl_8","ICNTL(8): scaling strategy (-2 to 7 or 77)","None",lu->id.ICNTL(8),&lu->id.ICNTL(8),PETSC_NULL);CHKERRQ(ierr);
738dcd589f8SShri Abhyankar   ierr = PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",lu->id.ICNTL(10),&lu->id.ICNTL(10),PETSC_NULL);CHKERRQ(ierr);
739dcd589f8SShri Abhyankar   ierr = PetscOptionsInt("-mat_mumps_icntl_11","ICNTL(11): statistics related to the linear system solved (via -ksp_view)","None",lu->id.ICNTL(11),&lu->id.ICNTL(11),PETSC_NULL);CHKERRQ(ierr);
740dcd589f8SShri Abhyankar   ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control: defines the ordering strategy with scaling constraints (0 to 3","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr);
741dcd589f8SShri Abhyankar   ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control: with or without ScaLAPACK","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr);
742dcd589f8SShri Abhyankar   ierr = PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage of estimated workspace increase","None",lu->id.ICNTL(14),&lu->id.ICNTL(14),PETSC_NULL);CHKERRQ(ierr);
743dcd589f8SShri Abhyankar   ierr = PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): Schur complement","None",lu->id.ICNTL(19),&lu->id.ICNTL(19),PETSC_NULL);CHKERRQ(ierr);
744dcd589f8SShri Abhyankar 
745dcd589f8SShri Abhyankar   ierr = PetscOptionsInt("-mat_mumps_icntl_22","ICNTL(22): in-core/out-of-core facility (0 or 1)","None",lu->id.ICNTL(22),&lu->id.ICNTL(22),PETSC_NULL);CHKERRQ(ierr);
746dcd589f8SShri Abhyankar   ierr = PetscOptionsInt("-mat_mumps_icntl_23","ICNTL(23): max size of the working memory (MB) that can allocate per processor","None",lu->id.ICNTL(23),&lu->id.ICNTL(23),PETSC_NULL);CHKERRQ(ierr);
747dcd589f8SShri Abhyankar   ierr = PetscOptionsInt("-mat_mumps_icntl_24","ICNTL(24): detection of null pivot rows (0 or 1)","None",lu->id.ICNTL(24),&lu->id.ICNTL(24),PETSC_NULL);CHKERRQ(ierr);
748dcd589f8SShri Abhyankar   ierr = PetscOptionsInt("-mat_mumps_icntl_25","ICNTL(25): computation of a null space basis","None",lu->id.ICNTL(25),&lu->id.ICNTL(25),PETSC_NULL);CHKERRQ(ierr);
749dcd589f8SShri Abhyankar   ierr = PetscOptionsInt("-mat_mumps_icntl_26","ICNTL(26): Schur options for right-hand side or solution vector","None",lu->id.ICNTL(26),&lu->id.ICNTL(26),PETSC_NULL);CHKERRQ(ierr);
750dcd589f8SShri Abhyankar   ierr = PetscOptionsInt("-mat_mumps_icntl_27","ICNTL(27): experimental parameter","None",lu->id.ICNTL(27),&lu->id.ICNTL(27),PETSC_NULL);CHKERRQ(ierr);
751292fb18eSBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_28","ICNTL(28): use 1 sequential (ictnl(7)) or 2 parallel (ictnl(29)) ordering","None",lu->id.ICNTL(28),&lu->id.ICNTL(28),PETSC_NULL);CHKERRQ(ierr);
752292fb18eSBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_29","ICNTL(29): parallel ordering 1 = ptscotch 2 = parmetis","None",lu->id.ICNTL(29),&lu->id.ICNTL(29),PETSC_NULL);CHKERRQ(ierr);
753dcd589f8SShri Abhyankar 
754dcd589f8SShri Abhyankar   ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",lu->id.CNTL(1),&lu->id.CNTL(1),PETSC_NULL);CHKERRQ(ierr);
755dcd589f8SShri Abhyankar   ierr = PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",lu->id.CNTL(2),&lu->id.CNTL(2),PETSC_NULL);CHKERRQ(ierr);
756dcd589f8SShri Abhyankar   ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",lu->id.CNTL(3),&lu->id.CNTL(3),PETSC_NULL);CHKERRQ(ierr);
757dcd589f8SShri Abhyankar   ierr = PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",lu->id.CNTL(4),&lu->id.CNTL(4),PETSC_NULL);CHKERRQ(ierr);
758dcd589f8SShri Abhyankar   ierr = PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",lu->id.CNTL(5),&lu->id.CNTL(5),PETSC_NULL);CHKERRQ(ierr);
759dcd589f8SShri Abhyankar   PetscOptionsEnd();
760dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
761dcd589f8SShri Abhyankar }
762dcd589f8SShri Abhyankar 
763dcd589f8SShri Abhyankar #undef __FUNCT__
764dcd589f8SShri Abhyankar #define __FUNCT__ "PetscInitializeMUMPS"
765f697e70eSHong Zhang PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS* mumps)
766dcd589f8SShri Abhyankar {
767dcd589f8SShri Abhyankar   PetscErrorCode  ierr;
768dcd589f8SShri Abhyankar 
769dcd589f8SShri Abhyankar   PetscFunctionBegin;
770f697e70eSHong Zhang   ierr = MPI_Comm_rank(((PetscObject)A)->comm, &mumps->myid);
771f697e70eSHong Zhang   ierr = MPI_Comm_size(((PetscObject)A)->comm,&mumps->size);CHKERRQ(ierr);
772f697e70eSHong Zhang   ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(mumps->comm_mumps));CHKERRQ(ierr);
773f697e70eSHong Zhang   mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps);
774f697e70eSHong Zhang 
775f697e70eSHong Zhang   mumps->id.job = JOB_INIT;
776f697e70eSHong Zhang   mumps->id.par = 1;  /* host participates factorizaton and solve */
777f697e70eSHong Zhang   mumps->id.sym = mumps->sym;
778dcd589f8SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
779f697e70eSHong Zhang   zmumps_c(&mumps->id);
780dcd589f8SShri Abhyankar #else
781f697e70eSHong Zhang   dmumps_c(&mumps->id);
782dcd589f8SShri Abhyankar #endif
783f697e70eSHong Zhang 
784f697e70eSHong Zhang   mumps->CleanUpMUMPS = PETSC_FALSE;
785f697e70eSHong Zhang   mumps->scat_rhs     = PETSC_NULL;
786f697e70eSHong Zhang   mumps->scat_sol     = PETSC_NULL;
787f697e70eSHong Zhang   mumps->nSolve       = 0;
788dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
789dcd589f8SShri Abhyankar }
790dcd589f8SShri Abhyankar 
791397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */
792397b6df1SKris Buschelman #undef __FUNCT__
793f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
7940481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
795b24902e0SBarry Smith {
796719d5645SBarry Smith   Mat_MUMPS          *lu = (Mat_MUMPS*)F->spptr;
797dcd589f8SShri Abhyankar   PetscErrorCode     ierr;
798bccb9932SShri Abhyankar   MatReuse           reuse;
79967877ebaSShri Abhyankar   Vec                b;
80067877ebaSShri Abhyankar   IS                 is_iden;
80167877ebaSShri Abhyankar   const PetscInt     M = A->rmap->N;
802397b6df1SKris Buschelman 
803397b6df1SKris Buschelman   PetscFunctionBegin;
804b24902e0SBarry Smith   lu->matstruc = DIFFERENT_NONZERO_PATTERN;
805dcd589f8SShri Abhyankar 
806dcd589f8SShri Abhyankar   /* Set MUMPS options */
807dcd589f8SShri Abhyankar   ierr = PetscSetMUMPSOptions(F,A);CHKERRQ(ierr);
808dcd589f8SShri Abhyankar 
809bccb9932SShri Abhyankar   reuse = MAT_INITIAL_MATRIX;
810bccb9932SShri Abhyankar   ierr = (*lu->ConvertToTriples)(A, 1, reuse, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
811dcd589f8SShri Abhyankar 
81267877ebaSShri Abhyankar   /* analysis phase */
81367877ebaSShri Abhyankar   /*----------------*/
8143d472b54SHong Zhang   lu->id.job = JOB_FACTSYMBOLIC;
81567877ebaSShri Abhyankar   lu->id.n = M;
81667877ebaSShri Abhyankar   switch (lu->id.ICNTL(18)){
81767877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
81867877ebaSShri Abhyankar     if (!lu->myid) {
81967877ebaSShri Abhyankar       lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
82067877ebaSShri Abhyankar       if (lu->id.ICNTL(6)>1){
82167877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
82267877ebaSShri Abhyankar         lu->id.a = (mumps_double_complex*)lu->val;
82367877ebaSShri Abhyankar #else
82467877ebaSShri Abhyankar         lu->id.a = lu->val;
82567877ebaSShri Abhyankar #endif
82667877ebaSShri Abhyankar       }
82767877ebaSShri Abhyankar     }
82867877ebaSShri Abhyankar     break;
82967877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
83067877ebaSShri Abhyankar     lu->id.nz_loc = lu->nz;
83167877ebaSShri Abhyankar     lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
83267877ebaSShri Abhyankar     if (lu->id.ICNTL(6)>1) {
83367877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
83467877ebaSShri Abhyankar       lu->id.a_loc = (mumps_double_complex*)lu->val;
83567877ebaSShri Abhyankar #else
83667877ebaSShri Abhyankar       lu->id.a_loc = lu->val;
83767877ebaSShri Abhyankar #endif
83867877ebaSShri Abhyankar     }
83967877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
84067877ebaSShri Abhyankar     if (!lu->myid){
84167877ebaSShri Abhyankar       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr);
84267877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
84367877ebaSShri Abhyankar     } else {
84467877ebaSShri Abhyankar       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr);
84567877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
84667877ebaSShri Abhyankar     }
84767877ebaSShri Abhyankar     ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr);
84867877ebaSShri Abhyankar     ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
84967877ebaSShri Abhyankar     ierr = VecSetFromOptions(b);CHKERRQ(ierr);
85067877ebaSShri Abhyankar 
85167877ebaSShri Abhyankar     ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr);
85267877ebaSShri Abhyankar     ierr = ISDestroy(is_iden);CHKERRQ(ierr);
85367877ebaSShri Abhyankar     ierr = VecDestroy(b);CHKERRQ(ierr);
85467877ebaSShri Abhyankar     break;
85567877ebaSShri Abhyankar     }
85667877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
85767877ebaSShri Abhyankar   zmumps_c(&lu->id);
85867877ebaSShri Abhyankar #else
85967877ebaSShri Abhyankar   dmumps_c(&lu->id);
86067877ebaSShri Abhyankar #endif
86167877ebaSShri 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));
86267877ebaSShri Abhyankar 
863719d5645SBarry Smith   F->ops->lufactornumeric  = MatFactorNumeric_MUMPS;
864dcd589f8SShri Abhyankar   F->ops->solve            = MatSolve_MUMPS;
865*51d5961aSHong Zhang   F->ops->solvetranspose   = MatSolveTranspose_MUMPS;
866b24902e0SBarry Smith   PetscFunctionReturn(0);
867b24902e0SBarry Smith }
868b24902e0SBarry Smith 
869450b117fSShri Abhyankar /* Note the Petsc r and c permutations are ignored */
870450b117fSShri Abhyankar #undef __FUNCT__
871450b117fSShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_BAIJMUMPS"
872450b117fSShri Abhyankar PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
873450b117fSShri Abhyankar {
874dcd589f8SShri Abhyankar 
875450b117fSShri Abhyankar   Mat_MUMPS       *lu = (Mat_MUMPS*)F->spptr;
876dcd589f8SShri Abhyankar   PetscErrorCode  ierr;
877bccb9932SShri Abhyankar   MatReuse        reuse;
87867877ebaSShri Abhyankar   Vec             b;
87967877ebaSShri Abhyankar   IS              is_iden;
88067877ebaSShri Abhyankar   const PetscInt  M = A->rmap->N;
881450b117fSShri Abhyankar 
882450b117fSShri Abhyankar   PetscFunctionBegin;
883450b117fSShri Abhyankar   lu->matstruc = DIFFERENT_NONZERO_PATTERN;
884dcd589f8SShri Abhyankar 
885dcd589f8SShri Abhyankar   /* Set MUMPS options */
886dcd589f8SShri Abhyankar   ierr = PetscSetMUMPSOptions(F,A);CHKERRQ(ierr);
887dcd589f8SShri Abhyankar 
888bccb9932SShri Abhyankar   reuse = MAT_INITIAL_MATRIX;
889bccb9932SShri Abhyankar   ierr = (*lu->ConvertToTriples)(A, 1, reuse, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
89067877ebaSShri Abhyankar 
89167877ebaSShri Abhyankar   /* analysis phase */
89267877ebaSShri Abhyankar   /*----------------*/
8933d472b54SHong Zhang   lu->id.job = JOB_FACTSYMBOLIC;
89467877ebaSShri Abhyankar   lu->id.n = M;
89567877ebaSShri Abhyankar   switch (lu->id.ICNTL(18)){
89667877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
89767877ebaSShri Abhyankar     if (!lu->myid) {
89867877ebaSShri Abhyankar       lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
89967877ebaSShri Abhyankar       if (lu->id.ICNTL(6)>1){
90067877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
90167877ebaSShri Abhyankar         lu->id.a = (mumps_double_complex*)lu->val;
90267877ebaSShri Abhyankar #else
90367877ebaSShri Abhyankar         lu->id.a = lu->val;
90467877ebaSShri Abhyankar #endif
90567877ebaSShri Abhyankar       }
90667877ebaSShri Abhyankar     }
90767877ebaSShri Abhyankar     break;
90867877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
90967877ebaSShri Abhyankar     lu->id.nz_loc = lu->nz;
91067877ebaSShri Abhyankar     lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
91167877ebaSShri Abhyankar     if (lu->id.ICNTL(6)>1) {
91267877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
91367877ebaSShri Abhyankar       lu->id.a_loc = (mumps_double_complex*)lu->val;
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);
93167877ebaSShri Abhyankar     ierr = ISDestroy(is_iden);CHKERRQ(ierr);
93267877ebaSShri Abhyankar     ierr = VecDestroy(b);CHKERRQ(ierr);
93367877ebaSShri Abhyankar     break;
93467877ebaSShri Abhyankar     }
93567877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
93667877ebaSShri Abhyankar   zmumps_c(&lu->id);
93767877ebaSShri Abhyankar #else
93867877ebaSShri Abhyankar   dmumps_c(&lu->id);
93967877ebaSShri Abhyankar #endif
94067877ebaSShri 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));
94167877ebaSShri Abhyankar 
942450b117fSShri Abhyankar   F->ops->lufactornumeric  = MatFactorNumeric_MUMPS;
943dcd589f8SShri Abhyankar   F->ops->solve            = MatSolve_MUMPS;
944*51d5961aSHong Zhang   F->ops->solvetranspose   = MatSolveTranspose_MUMPS;
945450b117fSShri Abhyankar   PetscFunctionReturn(0);
946450b117fSShri Abhyankar }
947b24902e0SBarry Smith 
948397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */
949397b6df1SKris Buschelman #undef __FUNCT__
95067877ebaSShri Abhyankar #define __FUNCT__ "MatCholeskyFactorSymbolic_MUMPS"
95167877ebaSShri Abhyankar PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
952b24902e0SBarry Smith {
9532792810eSHong Zhang   Mat_MUMPS          *lu = (Mat_MUMPS*)F->spptr;
954dcd589f8SShri Abhyankar   PetscErrorCode     ierr;
955bccb9932SShri Abhyankar   MatReuse           reuse;
95667877ebaSShri Abhyankar   Vec                b;
95767877ebaSShri Abhyankar   IS                 is_iden;
95867877ebaSShri Abhyankar   const PetscInt     M = A->rmap->N;
959397b6df1SKris Buschelman 
960397b6df1SKris Buschelman   PetscFunctionBegin;
961b24902e0SBarry Smith   lu->matstruc = DIFFERENT_NONZERO_PATTERN;
962dcd589f8SShri Abhyankar 
963dcd589f8SShri Abhyankar   /* Set MUMPS options */
964dcd589f8SShri Abhyankar   ierr = PetscSetMUMPSOptions(F,A);CHKERRQ(ierr);
965dcd589f8SShri Abhyankar 
966bccb9932SShri Abhyankar   reuse = MAT_INITIAL_MATRIX;
967bccb9932SShri Abhyankar   ierr = (*lu->ConvertToTriples)(A, 1 , reuse, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
968dcd589f8SShri Abhyankar 
96967877ebaSShri Abhyankar   /* analysis phase */
97067877ebaSShri Abhyankar   /*----------------*/
9713d472b54SHong Zhang   lu->id.job = JOB_FACTSYMBOLIC;
97267877ebaSShri Abhyankar   lu->id.n = M;
97367877ebaSShri Abhyankar   switch (lu->id.ICNTL(18)){
97467877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
97567877ebaSShri Abhyankar     if (!lu->myid) {
97667877ebaSShri Abhyankar       lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
97767877ebaSShri Abhyankar       if (lu->id.ICNTL(6)>1){
97867877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
97967877ebaSShri Abhyankar         lu->id.a = (mumps_double_complex*)lu->val;
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)
99167877ebaSShri Abhyankar       lu->id.a_loc = (mumps_double_complex*)lu->val;
99267877ebaSShri Abhyankar #else
99367877ebaSShri Abhyankar       lu->id.a_loc = lu->val;
99467877ebaSShri Abhyankar #endif
99567877ebaSShri Abhyankar     }
99667877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
99767877ebaSShri Abhyankar     if (!lu->myid){
99867877ebaSShri Abhyankar       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr);
99967877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
100067877ebaSShri Abhyankar     } else {
100167877ebaSShri Abhyankar       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr);
100267877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
100367877ebaSShri Abhyankar     }
100467877ebaSShri Abhyankar     ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr);
100567877ebaSShri Abhyankar     ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
100667877ebaSShri Abhyankar     ierr = VecSetFromOptions(b);CHKERRQ(ierr);
100767877ebaSShri Abhyankar 
100867877ebaSShri Abhyankar     ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr);
100967877ebaSShri Abhyankar     ierr = ISDestroy(is_iden);CHKERRQ(ierr);
101067877ebaSShri Abhyankar     ierr = VecDestroy(b);CHKERRQ(ierr);
101167877ebaSShri Abhyankar     break;
101267877ebaSShri Abhyankar     }
101367877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
101467877ebaSShri Abhyankar   zmumps_c(&lu->id);
101567877ebaSShri Abhyankar #else
101667877ebaSShri Abhyankar   dmumps_c(&lu->id);
101767877ebaSShri Abhyankar #endif
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 
10202792810eSHong Zhang   F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
1021dcd589f8SShri Abhyankar   F->ops->solve                 = MatSolve_MUMPS;
1022*51d5961aSHong Zhang   F->ops->solvetranspose        = MatSolve_MUMPS;
1023db4efbfdSBarry Smith #if !defined(PETSC_USE_COMPLEX)
1024dcd589f8SShri Abhyankar   (F)->ops->getinertia          = MatGetInertia_SBAIJMUMPS;
1025db4efbfdSBarry Smith #endif
1026b24902e0SBarry Smith   PetscFunctionReturn(0);
1027b24902e0SBarry Smith }
1028b24902e0SBarry Smith 
1029397b6df1SKris Buschelman #undef __FUNCT__
103064e6c443SBarry Smith #define __FUNCT__ "MatView_MUMPS"
103164e6c443SBarry Smith PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
103274ed9c26SBarry Smith {
1033f6c57405SHong Zhang   PetscErrorCode    ierr;
103464e6c443SBarry Smith   PetscBool         iascii;
103564e6c443SBarry Smith   PetscViewerFormat format;
103664e6c443SBarry Smith   Mat_MUMPS         *lu=(Mat_MUMPS*)A->spptr;
1037f6c57405SHong Zhang 
1038f6c57405SHong Zhang   PetscFunctionBegin;
103964e6c443SBarry Smith   /* check if matrix is mumps type */
104064e6c443SBarry Smith   if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0);
104164e6c443SBarry Smith 
104264e6c443SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
104364e6c443SBarry Smith   if (iascii) {
104464e6c443SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
104564e6c443SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO){
104664e6c443SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
104764e6c443SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                   %d \n",lu->id.sym);CHKERRQ(ierr);
104864e6c443SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):            %d \n",lu->id.par);CHKERRQ(ierr);
1049f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):         %d \n",lu->id.ICNTL(1));CHKERRQ(ierr);
1050f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg): %d \n",lu->id.ICNTL(2));CHKERRQ(ierr);
1051f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):   %d \n",lu->id.ICNTL(3));CHKERRQ(ierr);
1052f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):        %d \n",lu->id.ICNTL(4));CHKERRQ(ierr);
1053f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):         %d \n",lu->id.ICNTL(5));CHKERRQ(ierr);
1054f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):        %d \n",lu->id.ICNTL(6));CHKERRQ(ierr);
1055d06efebcSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (sequentia matrix ordering):%d \n",lu->id.ICNTL(7));CHKERRQ(ierr);
1056f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scalling strategy):        %d \n",lu->id.ICNTL(8));CHKERRQ(ierr);
1057f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements):  %d \n",lu->id.ICNTL(10));CHKERRQ(ierr);
1058f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):          %d \n",lu->id.ICNTL(11));CHKERRQ(ierr);
105934ed7027SBarry Smith       if (lu->id.ICNTL(11)>0) {
106034ed7027SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(4) (inf norm of input mat):        %g\n",lu->id.RINFOG(4));CHKERRQ(ierr);
106134ed7027SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(5) (inf norm of solution):         %g\n",lu->id.RINFOG(5));CHKERRQ(ierr);
106234ed7027SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(6) (inf norm of residual):         %g\n",lu->id.RINFOG(6));CHKERRQ(ierr);
106334ed7027SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",lu->id.RINFOG(7),lu->id.RINFOG(8));CHKERRQ(ierr);
106434ed7027SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(9) (error estimate):               %g \n",lu->id.RINFOG(9));CHKERRQ(ierr);
106534ed7027SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",lu->id.RINFOG(10),lu->id.RINFOG(11));CHKERRQ(ierr);
1066f6c57405SHong Zhang       }
1067f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",lu->id.ICNTL(12));CHKERRQ(ierr);
1068f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",lu->id.ICNTL(13));CHKERRQ(ierr);
1069f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr);
1070f6c57405SHong Zhang       /* ICNTL(15-17) not used */
1071f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",lu->id.ICNTL(18));CHKERRQ(ierr);
1072f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Shur complement info):                       %d \n",lu->id.ICNTL(19));CHKERRQ(ierr);
1073f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",lu->id.ICNTL(20));CHKERRQ(ierr);
1074f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d \n",lu->id.ICNTL(21));CHKERRQ(ierr);
1075c0165424SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",lu->id.ICNTL(22));CHKERRQ(ierr);
1076c0165424SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",lu->id.ICNTL(23));CHKERRQ(ierr);
1077c0165424SHong Zhang 
1078c0165424SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",lu->id.ICNTL(24));CHKERRQ(ierr);
1079c0165424SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",lu->id.ICNTL(25));CHKERRQ(ierr);
1080c0165424SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",lu->id.ICNTL(26));CHKERRQ(ierr);
1081c0165424SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",lu->id.ICNTL(27));CHKERRQ(ierr);
1082d06efebcSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(28) (Use parallel or sequential ordering):        %d \n",lu->id.ICNTL(28));CHKERRQ(ierr);
1083d06efebcSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(29) (Parallel ordering):                          %d \n",lu->id.ICNTL(29));CHKERRQ(ierr);
1084f6c57405SHong Zhang 
1085f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",lu->id.CNTL(1));CHKERRQ(ierr);
1086f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr);
1087f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",lu->id.CNTL(3));CHKERRQ(ierr);
1088f6c57405SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",lu->id.CNTL(4));CHKERRQ(ierr);
1089c0165424SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",lu->id.CNTL(5));CHKERRQ(ierr);
1090f6c57405SHong Zhang 
1091f6c57405SHong Zhang       /* infomation local to each processor */
109234ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);
109334ed7027SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr);
109434ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
109534ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);
109634ed7027SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr);
109734ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
109834ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);
109934ed7027SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr);
110034ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1101f6c57405SHong Zhang 
110234ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr);
110334ed7027SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] %d \n",lu->myid,lu->id.INFO(15));CHKERRQ(ierr);
110434ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1105f6c57405SHong Zhang 
110634ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr);
110734ed7027SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",lu->myid,lu->id.INFO(16));CHKERRQ(ierr);
110834ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1109f6c57405SHong Zhang 
111034ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);
111134ed7027SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",lu->myid,lu->id.INFO(23));CHKERRQ(ierr);
111234ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1113f6c57405SHong Zhang 
1114f6c57405SHong Zhang       if (!lu->myid){ /* information from the host */
1115f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr);
1116f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr);
1117f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr);
1118f6c57405SHong Zhang 
1119f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr);
1120f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr);
1121f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr);
1122f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr);
1123f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr);
1124f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr);
1125f6c57405SHong 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);
1126f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",lu->id.INFOG(10));CHKERRQ(ierr);
1127f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",lu->id.INFOG(11));CHKERRQ(ierr);
1128f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr);
1129f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr);
1130f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr);
1131f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr);
1132f6c57405SHong 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);
1133f6c57405SHong 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);
1134f6c57405SHong 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);
1135f6c57405SHong 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);
1136f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr);
1137f6c57405SHong 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);
1138f6c57405SHong 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);
1139f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",lu->id.INFOG(23));CHKERRQ(ierr);
1140f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",lu->id.INFOG(24));CHKERRQ(ierr);
1141f6c57405SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",lu->id.INFOG(25));CHKERRQ(ierr);
1142f6c57405SHong Zhang       }
1143f6c57405SHong Zhang     }
1144cb828f0fSHong Zhang   }
1145f6c57405SHong Zhang   PetscFunctionReturn(0);
1146f6c57405SHong Zhang }
1147f6c57405SHong Zhang 
114835bd34faSBarry Smith #undef __FUNCT__
114935bd34faSBarry Smith #define __FUNCT__ "MatGetInfo_MUMPS"
115035bd34faSBarry Smith PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
115135bd34faSBarry Smith {
1152cb828f0fSHong Zhang   Mat_MUMPS  *mumps =(Mat_MUMPS*)A->spptr;
115335bd34faSBarry Smith 
115435bd34faSBarry Smith   PetscFunctionBegin;
115535bd34faSBarry Smith   info->block_size        = 1.0;
1156cb828f0fSHong Zhang   info->nz_allocated      = mumps->id.INFOG(20);
1157cb828f0fSHong Zhang   info->nz_used           = mumps->id.INFOG(20);
115835bd34faSBarry Smith   info->nz_unneeded       = 0.0;
115935bd34faSBarry Smith   info->assemblies        = 0.0;
116035bd34faSBarry Smith   info->mallocs           = 0.0;
116135bd34faSBarry Smith   info->memory            = 0.0;
116235bd34faSBarry Smith   info->fill_ratio_given  = 0;
116335bd34faSBarry Smith   info->fill_ratio_needed = 0;
116435bd34faSBarry Smith   info->factor_mallocs    = 0;
116535bd34faSBarry Smith   PetscFunctionReturn(0);
116635bd34faSBarry Smith }
116735bd34faSBarry Smith 
11685ccb76cbSHong Zhang /* -------------------------------------------------------------------------------------------*/
11695ccb76cbSHong Zhang #undef __FUNCT__
11705ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl_MUMPS"
11715ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival)
11725ccb76cbSHong Zhang {
11735ccb76cbSHong Zhang   Mat_MUMPS *lu =(Mat_MUMPS*)F->spptr;
11745ccb76cbSHong Zhang 
11755ccb76cbSHong Zhang   PetscFunctionBegin;
11765ccb76cbSHong Zhang   lu->id.ICNTL(icntl) = ival;
11775ccb76cbSHong Zhang   PetscFunctionReturn(0);
11785ccb76cbSHong Zhang }
11795ccb76cbSHong Zhang 
11805ccb76cbSHong Zhang #undef __FUNCT__
11815ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl"
11825ccb76cbSHong Zhang /*@
11835ccb76cbSHong Zhang   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()
11845ccb76cbSHong Zhang 
11855ccb76cbSHong Zhang    Logically Collective on Mat
11865ccb76cbSHong Zhang 
11875ccb76cbSHong Zhang    Input Parameters:
11885ccb76cbSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
11895ccb76cbSHong Zhang .  icntl - index of MUMPS parameter array ICNTL()
11905ccb76cbSHong Zhang -  ival - value of MUMPS ICNTL(icntl)
11915ccb76cbSHong Zhang 
11925ccb76cbSHong Zhang   Options Database:
11935ccb76cbSHong Zhang .   -mat_mumps_icntl_<icntl> <ival>
11945ccb76cbSHong Zhang 
11955ccb76cbSHong Zhang    Level: beginner
11965ccb76cbSHong Zhang 
11975ccb76cbSHong Zhang    References: MUMPS Users' Guide
11985ccb76cbSHong Zhang 
11995ccb76cbSHong Zhang .seealso: MatGetFactor()
12005ccb76cbSHong Zhang @*/
12015ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival)
12025ccb76cbSHong Zhang {
12035ccb76cbSHong Zhang   PetscErrorCode ierr;
12045ccb76cbSHong Zhang 
12055ccb76cbSHong Zhang   PetscFunctionBegin;
12065ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
12075ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,ival,3);
12085ccb76cbSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr);
12095ccb76cbSHong Zhang   PetscFunctionReturn(0);
12105ccb76cbSHong Zhang }
12115ccb76cbSHong Zhang 
121224b6179bSKris Buschelman /*MC
12132692d6eeSBarry Smith   MATSOLVERMUMPS -  A matrix type providing direct solvers (LU and Cholesky) for
121424b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
121524b6179bSKris Buschelman 
121641c8de11SBarry Smith   Works with MATAIJ and MATSBAIJ matrices
121724b6179bSKris Buschelman 
121824b6179bSKris Buschelman   Options Database Keys:
1219fb8376fbSHong Zhang + -mat_mumps_icntl_4 <0,...,4> - print level
122024b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
122164e6c443SBarry Smith . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guidec)
122224b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
122324b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
122494b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
122524b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
122624b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
122724b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
122824b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
122924b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
123024b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
123124b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
123224b6179bSKris Buschelman 
123324b6179bSKris Buschelman   Level: beginner
123424b6179bSKris Buschelman 
123541c8de11SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage
123641c8de11SBarry Smith 
123724b6179bSKris Buschelman M*/
123824b6179bSKris Buschelman 
12392877fffaSHong Zhang EXTERN_C_BEGIN
124035bd34faSBarry Smith #undef __FUNCT__
124135bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_mumps"
124235bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type)
124335bd34faSBarry Smith {
124435bd34faSBarry Smith   PetscFunctionBegin;
12452692d6eeSBarry Smith   *type = MATSOLVERMUMPS;
124635bd34faSBarry Smith   PetscFunctionReturn(0);
124735bd34faSBarry Smith }
124835bd34faSBarry Smith EXTERN_C_END
124935bd34faSBarry Smith 
125035bd34faSBarry Smith EXTERN_C_BEGIN
1251bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI AIJ matrices */
12522877fffaSHong Zhang #undef __FUNCT__
1253bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_aij_mumps"
1254bccb9932SShri Abhyankar PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
12552877fffaSHong Zhang {
12562877fffaSHong Zhang   Mat            B;
12572877fffaSHong Zhang   PetscErrorCode ierr;
12582877fffaSHong Zhang   Mat_MUMPS      *mumps;
1259ace3abfcSBarry Smith   PetscBool      isSeqAIJ;
12602877fffaSHong Zhang 
12612877fffaSHong Zhang   PetscFunctionBegin;
12622877fffaSHong Zhang   /* Create the factorization matrix */
1263bccb9932SShri Abhyankar   ierr = PetscTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr);
12642877fffaSHong Zhang   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
12652877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
12662877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1267bccb9932SShri Abhyankar   if (isSeqAIJ) {
12682877fffaSHong Zhang     ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
1269bccb9932SShri Abhyankar   } else {
1270bccb9932SShri Abhyankar     ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
1271bccb9932SShri Abhyankar   }
12722877fffaSHong Zhang 
127316ebf90aSShri Abhyankar   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
12742877fffaSHong Zhang   B->ops->view             = MatView_MUMPS;
127535bd34faSBarry Smith   B->ops->getinfo          = MatGetInfo_MUMPS;
127635bd34faSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
12775ccb76cbSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl_MUMPS",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
1278450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
1279450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
1280d5f3da31SBarry Smith     B->factortype = MAT_FACTOR_LU;
1281bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
1282bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
1283746480a1SHong Zhang     mumps->sym = 0;
1284dcd589f8SShri Abhyankar   } else {
128567877ebaSShri Abhyankar     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
1286450b117fSShri Abhyankar     B->factortype = MAT_FACTOR_CHOLESKY;
1287bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
1288bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
12896fdc2a6dSBarry Smith     if (A->spd_set && A->spd) mumps->sym = 1;
12906fdc2a6dSBarry Smith     else                      mumps->sym = 2;
1291450b117fSShri Abhyankar   }
12922877fffaSHong Zhang 
12932877fffaSHong Zhang   mumps->isAIJ        = PETSC_TRUE;
12942877fffaSHong Zhang   mumps->MatDestroy   = B->ops->destroy;
12952877fffaSHong Zhang   B->ops->destroy     = MatDestroy_MUMPS;
12962877fffaSHong Zhang   B->spptr            = (void*)mumps;
1297f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1298746480a1SHong Zhang 
12992877fffaSHong Zhang   *F = B;
13002877fffaSHong Zhang   PetscFunctionReturn(0);
13012877fffaSHong Zhang }
13022877fffaSHong Zhang EXTERN_C_END
13032877fffaSHong Zhang 
1304bccb9932SShri Abhyankar 
13052877fffaSHong Zhang EXTERN_C_BEGIN
1306bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI SBAIJ matrices */
13072877fffaSHong Zhang #undef __FUNCT__
1308bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_sbaij_mumps"
1309bccb9932SShri Abhyankar PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
13102877fffaSHong Zhang {
13112877fffaSHong Zhang   Mat            B;
13122877fffaSHong Zhang   PetscErrorCode ierr;
13132877fffaSHong Zhang   Mat_MUMPS      *mumps;
1314ace3abfcSBarry Smith   PetscBool      isSeqSBAIJ;
13152877fffaSHong Zhang 
13162877fffaSHong Zhang   PetscFunctionBegin;
1317e7e72b3dSBarry Smith   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
1318bccb9932SShri 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");
1319bccb9932SShri Abhyankar   ierr = PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
13202877fffaSHong Zhang   /* Create the factorization matrix */
13212877fffaSHong Zhang   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
13222877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
13232877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
132416ebf90aSShri Abhyankar   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
1325bccb9932SShri Abhyankar   if (isSeqSBAIJ) {
1326bccb9932SShri Abhyankar     ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
132716ebf90aSShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
1328dcd589f8SShri Abhyankar   } else {
1329bccb9932SShri Abhyankar     ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
1330bccb9932SShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
1331bccb9932SShri Abhyankar   }
1332bccb9932SShri Abhyankar 
133367877ebaSShri Abhyankar   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
1334bccb9932SShri Abhyankar   B->ops->view                   = MatView_MUMPS;
1335bccb9932SShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
1336f250808bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl",MatMumpsSetIcntl);CHKERRQ(ierr);
1337f4762488SHong Zhang   B->factortype                  = MAT_FACTOR_CHOLESKY;
13386fdc2a6dSBarry Smith   if (A->spd_set && A->spd) mumps->sym = 1;
13396fdc2a6dSBarry Smith   else                      mumps->sym = 2;
1340a214ac2aSShri Abhyankar 
1341bccb9932SShri Abhyankar   mumps->isAIJ        = PETSC_FALSE;
1342f3c0ef26SHong Zhang   mumps->MatDestroy   = B->ops->destroy;
1343f3c0ef26SHong Zhang   B->ops->destroy     = MatDestroy_MUMPS;
13442877fffaSHong Zhang   B->spptr            = (void*)mumps;
1345f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1346746480a1SHong Zhang 
13472877fffaSHong Zhang   *F = B;
13482877fffaSHong Zhang   PetscFunctionReturn(0);
13492877fffaSHong Zhang }
13502877fffaSHong Zhang EXTERN_C_END
135197969023SHong Zhang 
1352450b117fSShri Abhyankar EXTERN_C_BEGIN
1353450b117fSShri Abhyankar #undef __FUNCT__
1354bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_baij_mumps"
1355bccb9932SShri Abhyankar PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
135667877ebaSShri Abhyankar {
135767877ebaSShri Abhyankar   Mat            B;
135867877ebaSShri Abhyankar   PetscErrorCode ierr;
135967877ebaSShri Abhyankar   Mat_MUMPS      *mumps;
1360ace3abfcSBarry Smith   PetscBool      isSeqBAIJ;
136167877ebaSShri Abhyankar 
136267877ebaSShri Abhyankar   PetscFunctionBegin;
136367877ebaSShri Abhyankar   /* Create the factorization matrix */
1364bccb9932SShri Abhyankar   ierr = PetscTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr);
136567877ebaSShri Abhyankar   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
136667877ebaSShri Abhyankar   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
136767877ebaSShri Abhyankar   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1368bccb9932SShri Abhyankar   if (isSeqBAIJ) {
136967877ebaSShri Abhyankar     ierr = MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,PETSC_NULL);CHKERRQ(ierr);
1370bccb9932SShri Abhyankar   } else {
137167877ebaSShri Abhyankar     ierr = MatMPIBAIJSetPreallocation(B,A->rmap->bs,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
1372bccb9932SShri Abhyankar   }
1373450b117fSShri Abhyankar 
137467877ebaSShri Abhyankar   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
1375450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
1376450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
1377450b117fSShri Abhyankar     B->factortype = MAT_FACTOR_LU;
1378bccb9932SShri Abhyankar     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
1379bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
1380746480a1SHong Zhang     mumps->sym = 0;
1381746480a1SHong Zhang   } else {
1382746480a1SHong Zhang     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n");
1383450b117fSShri Abhyankar   }
1384bccb9932SShri Abhyankar 
1385450b117fSShri Abhyankar   B->ops->view             = MatView_MUMPS;
1386450b117fSShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
13875ccb76cbSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl_MUMPS",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
1388450b117fSShri Abhyankar 
1389450b117fSShri Abhyankar   mumps->isAIJ        = PETSC_TRUE;
1390450b117fSShri Abhyankar   mumps->MatDestroy   = B->ops->destroy;
1391450b117fSShri Abhyankar   B->ops->destroy     = MatDestroy_MUMPS;
1392450b117fSShri Abhyankar   B->spptr            = (void*)mumps;
1393f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1394746480a1SHong Zhang 
1395450b117fSShri Abhyankar   *F = B;
1396450b117fSShri Abhyankar   PetscFunctionReturn(0);
1397450b117fSShri Abhyankar }
1398450b117fSShri Abhyankar EXTERN_C_END
1399a214ac2aSShri Abhyankar 
1400