xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision 185f65962d6549214bfeb8510bb5df406dbc8f62)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
21c2a3de1SBarry Smith 
3397b6df1SKris Buschelman /*
4c2b5dc30SHong Zhang     Provides an interface to the MUMPS sparse solver
5397b6df1SKris Buschelman */
64e1dd20eSHong Zhang #include "../src/mat/impls/aij/seq/aij.h"  /*I  "petscmat.h"  I*/
77c4f633dSBarry Smith #include "../src/mat/impls/aij/mpi/mpiaij.h"
87c4f633dSBarry Smith #include "../src/mat/impls/sbaij/seq/sbaij.h"
97c4f633dSBarry Smith #include "../src/mat/impls/sbaij/mpi/mpisbaij.h"
10450b117fSShri Abhyankar #include "../src/mat/impls/baij/seq/baij.h"
11450b117fSShri Abhyankar #include "../src/mat/impls/baij/mpi/mpibaij.h"
12397b6df1SKris Buschelman 
13397b6df1SKris Buschelman EXTERN_C_BEGIN
14397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
15397b6df1SKris Buschelman #include "zmumps_c.h"
16397b6df1SKris Buschelman #else
17397b6df1SKris Buschelman #include "dmumps_c.h"
18397b6df1SKris Buschelman #endif
19397b6df1SKris Buschelman EXTERN_C_END
20397b6df1SKris Buschelman #define JOB_INIT -1
21397b6df1SKris Buschelman #define JOB_END -2
22397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */
23397b6df1SKris Buschelman #define ICNTL(I) icntl[(I)-1]
24397b6df1SKris Buschelman #define CNTL(I) cntl[(I)-1]
25397b6df1SKris Buschelman #define INFOG(I) infog[(I)-1]
26a7aca84bSHong Zhang #define INFO(I) info[(I)-1]
27397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1]
28adc1d99fSHong Zhang #define RINFO(I) rinfo[(I)-1]
29397b6df1SKris Buschelman 
30397b6df1SKris Buschelman typedef struct {
31397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
32397b6df1SKris Buschelman   ZMUMPS_STRUC_C id;
33397b6df1SKris Buschelman #else
34397b6df1SKris Buschelman   DMUMPS_STRUC_C id;
35397b6df1SKris Buschelman #endif
36397b6df1SKris Buschelman   MatStructure   matstruc;
37c1490034SHong Zhang   PetscMPIInt    myid,size;
3816ebf90aSShri Abhyankar   PetscInt       *irn,*jcn,nz,sym,nSolve;
39397b6df1SKris Buschelman   PetscScalar    *val;
40397b6df1SKris Buschelman   MPI_Comm       comm_mumps;
41329ec9b3SHong Zhang   VecScatter     scat_rhs, scat_sol;
42cb828f0fSHong Zhang   PetscTruth     isAIJ,CleanUpMUMPS,mumpsview;
43329ec9b3SHong Zhang   Vec            b_seq,x_seq;
4467334b25SHong Zhang   PetscErrorCode (*MatDestroy)(Mat);
45bccb9932SShri Abhyankar   PetscErrorCode (*ConvertToTriples)(Mat, int, MatReuse, int*, int**, int**, PetscScalar**);
46f0c56d0fSKris Buschelman } Mat_MUMPS;
47f0c56d0fSKris Buschelman 
48dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);
49b24902e0SBarry Smith 
5067877ebaSShri Abhyankar 
5167877ebaSShri Abhyankar /* MatConvertToTriples_A_B */
5267877ebaSShri Abhyankar /*convert Petsc matrix to triples: row[nz], col[nz], val[nz] */
53397b6df1SKris Buschelman /*
54397b6df1SKris Buschelman   input:
5567877ebaSShri Abhyankar     A       - matrix in aij,baij or sbaij (bs=1) format
56397b6df1SKris Buschelman     shift   - 0: C style output triple; 1: Fortran style output triple.
57bccb9932SShri Abhyankar     reuse   - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple
58bccb9932SShri Abhyankar               MAT_REUSE_MATRIX:   only the values in v array are updated
59397b6df1SKris Buschelman   output:
60397b6df1SKris Buschelman     nnz     - dim of r, c, and v (number of local nonzero entries of A)
61397b6df1SKris Buschelman     r, c, v - row and col index, matrix values (matrix triples)
62397b6df1SKris Buschelman  */
6316ebf90aSShri Abhyankar 
6416ebf90aSShri Abhyankar #undef __FUNCT__
6516ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqaij"
66bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
67b24902e0SBarry Smith {
68*185f6596SHong Zhang   const PetscInt   *ai,*aj,*ajj,M=A->rmap->n;
6967877ebaSShri Abhyankar   PetscInt         nz,rnz,i,j;
70dfbe8321SBarry Smith   PetscErrorCode   ierr;
71c1490034SHong Zhang   PetscInt         *row,*col;
7216ebf90aSShri Abhyankar   Mat_SeqAIJ       *aa=(Mat_SeqAIJ*)A->data;
73397b6df1SKris Buschelman 
74397b6df1SKris Buschelman   PetscFunctionBegin;
7516ebf90aSShri Abhyankar   *v=aa->a;
76bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX){
7716ebf90aSShri Abhyankar     nz = aa->nz; ai = aa->i; aj = aa->j;
7816ebf90aSShri Abhyankar     *nnz = nz;
79*185f6596SHong Zhang     ierr = PetscMalloc(2*nz*sizeof(PetscInt), &row);CHKERRQ(ierr);
80*185f6596SHong Zhang     col  = row + nz;
81*185f6596SHong Zhang 
8216ebf90aSShri Abhyankar     nz = 0;
8316ebf90aSShri Abhyankar     for(i=0; i<M; i++) {
8416ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
8567877ebaSShri Abhyankar       ajj = aj + ai[i];
8667877ebaSShri Abhyankar       for(j=0; j<rnz; j++) {
8767877ebaSShri Abhyankar 	row[nz] = i+shift; col[nz++] = ajj[j] + shift;
8816ebf90aSShri Abhyankar       }
8916ebf90aSShri Abhyankar     }
9016ebf90aSShri Abhyankar     *r = row; *c = col;
9116ebf90aSShri Abhyankar   }
9216ebf90aSShri Abhyankar   PetscFunctionReturn(0);
9316ebf90aSShri Abhyankar }
94397b6df1SKris Buschelman 
9516ebf90aSShri Abhyankar #undef __FUNCT__
9667877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqbaij_seqaij"
97bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
9867877ebaSShri Abhyankar {
9967877ebaSShri Abhyankar   Mat_SeqBAIJ        *aa=(Mat_SeqBAIJ*)A->data;
10067877ebaSShri Abhyankar   const PetscInt     *ai,*aj,*ajj,bs=A->rmap->bs,bs2=aa->bs2,M=A->rmap->N/bs;
10167877ebaSShri Abhyankar   PetscInt           nz,idx=0,rnz,i,j,k,m,ii;
10267877ebaSShri Abhyankar   PetscErrorCode     ierr;
10367877ebaSShri Abhyankar   PetscInt           *row,*col;
10467877ebaSShri Abhyankar 
10567877ebaSShri Abhyankar   PetscFunctionBegin;
106cf3759fdSShri Abhyankar   *v = aa->a;
107bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX){
108cf3759fdSShri Abhyankar     ai = aa->i; aj = aa->j;
10967877ebaSShri Abhyankar     nz = bs2*aa->nz;
11067877ebaSShri Abhyankar     *nnz = nz;
111*185f6596SHong Zhang     ierr = PetscMalloc(2*nz*sizeof(PetscInt), &row);CHKERRQ(ierr);
112*185f6596SHong Zhang     col  = row + nz;
113*185f6596SHong Zhang 
11467877ebaSShri Abhyankar     for(i=0; i<M; i++) {
11567877ebaSShri Abhyankar       ii = 0;
11667877ebaSShri Abhyankar       ajj = aj + ai[i];
11767877ebaSShri Abhyankar       rnz = ai[i+1] - ai[i];
11867877ebaSShri Abhyankar       for(k=0; k<rnz; k++) {
11967877ebaSShri Abhyankar 	for(j=0; j<bs; j++) {
12067877ebaSShri Abhyankar 	  for(m=0; m<bs; m++) {
12167877ebaSShri Abhyankar 	    row[idx]     = i*bs + m + shift;
122cf3759fdSShri Abhyankar 	    col[idx++]   = bs*(ajj[k]) + j + shift;
12367877ebaSShri Abhyankar 	  }
12467877ebaSShri Abhyankar 	}
12567877ebaSShri Abhyankar       }
12667877ebaSShri Abhyankar     }
127cf3759fdSShri Abhyankar     *r = row; *c = col;
12867877ebaSShri Abhyankar   }
12967877ebaSShri Abhyankar   PetscFunctionReturn(0);
13067877ebaSShri Abhyankar }
13167877ebaSShri Abhyankar 
13267877ebaSShri Abhyankar #undef __FUNCT__
13316ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqsbaij_seqsbaij"
134bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
13516ebf90aSShri Abhyankar {
13667877ebaSShri Abhyankar   const PetscInt   *ai, *aj,*ajj,M=A->rmap->n;
13767877ebaSShri Abhyankar   PetscInt         nz,rnz,i,j;
13816ebf90aSShri Abhyankar   PetscErrorCode   ierr;
13916ebf90aSShri Abhyankar   PetscInt         *row,*col;
14016ebf90aSShri Abhyankar   Mat_SeqSBAIJ     *aa=(Mat_SeqSBAIJ*)A->data;
14116ebf90aSShri Abhyankar 
14216ebf90aSShri Abhyankar   PetscFunctionBegin;
143bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX){
14416ebf90aSShri Abhyankar     nz = aa->nz;ai=aa->i; aj=aa->j;*v=aa->a;
14516ebf90aSShri Abhyankar     *nnz = nz;
146*185f6596SHong Zhang     ierr = PetscMalloc(2*nz*sizeof(PetscInt), &row);CHKERRQ(ierr);
147*185f6596SHong Zhang     col  = row + nz;
148*185f6596SHong Zhang 
14916ebf90aSShri Abhyankar     nz = 0;
15016ebf90aSShri Abhyankar     for(i=0; i<M; i++) {
15116ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
15267877ebaSShri Abhyankar       ajj = aj + ai[i];
15367877ebaSShri Abhyankar       for(j=0; j<rnz; j++) {
15467877ebaSShri Abhyankar 	row[nz] = i+shift; col[nz++] = ajj[j] + shift;
15516ebf90aSShri Abhyankar       }
15616ebf90aSShri Abhyankar     }
15716ebf90aSShri Abhyankar     *r = row; *c = col;
15816ebf90aSShri Abhyankar   }
15916ebf90aSShri Abhyankar   PetscFunctionReturn(0);
16016ebf90aSShri Abhyankar }
16116ebf90aSShri Abhyankar 
16216ebf90aSShri Abhyankar #undef __FUNCT__
16316ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqsbaij"
164bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
16516ebf90aSShri Abhyankar {
16667877ebaSShri Abhyankar   const PetscInt     *ai,*aj,*ajj,*adiag,M=A->rmap->n;
16767877ebaSShri Abhyankar   PetscInt           nz,rnz,i,j;
16867877ebaSShri Abhyankar   const PetscScalar  *av,*v1;
16916ebf90aSShri Abhyankar   PetscScalar        *val;
17016ebf90aSShri Abhyankar   PetscErrorCode     ierr;
17116ebf90aSShri Abhyankar   PetscInt           *row,*col;
17216ebf90aSShri Abhyankar   Mat_SeqSBAIJ       *aa=(Mat_SeqSBAIJ*)A->data;
17316ebf90aSShri Abhyankar 
17416ebf90aSShri Abhyankar   PetscFunctionBegin;
17516ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j;av=aa->a;
17616ebf90aSShri Abhyankar   adiag=aa->diag;
177bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX){
17816ebf90aSShri Abhyankar     nz = M + (aa->nz-M)/2;
17916ebf90aSShri Abhyankar     *nnz = nz;
180*185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
181*185f6596SHong Zhang     col  = row + nz;
182*185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
183*185f6596SHong Zhang 
18416ebf90aSShri Abhyankar     nz = 0;
18516ebf90aSShri Abhyankar     for(i=0; i<M; i++) {
18616ebf90aSShri Abhyankar       rnz = ai[i+1] - adiag[i];
18767877ebaSShri Abhyankar       ajj  = aj + adiag[i];
188cf3759fdSShri Abhyankar       v1   = av + adiag[i];
18967877ebaSShri Abhyankar       for(j=0; j<rnz; j++) {
19067877ebaSShri Abhyankar 	row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
19116ebf90aSShri Abhyankar       }
19216ebf90aSShri Abhyankar     }
19316ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
194397b6df1SKris Buschelman   } else {
19516ebf90aSShri Abhyankar     nz = 0; val = *v;
19616ebf90aSShri Abhyankar     for(i=0; i <M; i++) {
19716ebf90aSShri Abhyankar       rnz = ai[i+1] - adiag[i];
19867877ebaSShri Abhyankar       ajj = aj + adiag[i];
19967877ebaSShri Abhyankar       v1  = av + adiag[i];
20067877ebaSShri Abhyankar       for(j=0; j<rnz; j++) {
20167877ebaSShri Abhyankar 	val[nz++] = v1[j];
20216ebf90aSShri Abhyankar       }
20316ebf90aSShri Abhyankar     }
20416ebf90aSShri Abhyankar   }
20516ebf90aSShri Abhyankar   PetscFunctionReturn(0);
20616ebf90aSShri Abhyankar }
20716ebf90aSShri Abhyankar 
20816ebf90aSShri Abhyankar #undef __FUNCT__
20916ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpisbaij_mpisbaij"
210bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
21116ebf90aSShri Abhyankar {
21216ebf90aSShri Abhyankar   const PetscInt     *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
21316ebf90aSShri Abhyankar   PetscErrorCode     ierr;
21416ebf90aSShri Abhyankar   PetscInt           rstart,nz,i,j,jj,irow,countA,countB;
21516ebf90aSShri Abhyankar   PetscInt           *row,*col;
21616ebf90aSShri Abhyankar   const PetscScalar  *av, *bv,*v1,*v2;
21716ebf90aSShri Abhyankar   PetscScalar        *val;
218397b6df1SKris Buschelman   Mat_MPISBAIJ       *mat =  (Mat_MPISBAIJ*)A->data;
219397b6df1SKris Buschelman   Mat_SeqSBAIJ       *aa=(Mat_SeqSBAIJ*)(mat->A)->data;
220397b6df1SKris Buschelman   Mat_SeqBAIJ        *bb=(Mat_SeqBAIJ*)(mat->B)->data;
22116ebf90aSShri Abhyankar 
22216ebf90aSShri Abhyankar   PetscFunctionBegin;
223d0f46423SBarry Smith   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
224397b6df1SKris Buschelman   garray = mat->garray;
225397b6df1SKris Buschelman   av=aa->a; bv=bb->a;
226397b6df1SKris Buschelman 
227bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX){
22816ebf90aSShri Abhyankar     nz = aa->nz + bb->nz;
22916ebf90aSShri Abhyankar     *nnz = nz;
230*185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
231*185f6596SHong Zhang     col  = row + nz;
232*185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
233*185f6596SHong Zhang 
234397b6df1SKris Buschelman     *r = row; *c = col; *v = val;
235397b6df1SKris Buschelman   } else {
236397b6df1SKris Buschelman     row = *r; col = *c; val = *v;
237397b6df1SKris Buschelman   }
238397b6df1SKris Buschelman 
239028e57e8SHong Zhang   jj = 0; irow = rstart;
240397b6df1SKris Buschelman   for ( i=0; i<m; i++ ) {
241397b6df1SKris Buschelman     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
242397b6df1SKris Buschelman     countA = ai[i+1] - ai[i];
243397b6df1SKris Buschelman     countB = bi[i+1] - bi[i];
244397b6df1SKris Buschelman     bjj    = bj + bi[i];
24516ebf90aSShri Abhyankar     v1     = av + ai[i];
24616ebf90aSShri Abhyankar     v2     = bv + bi[i];
247397b6df1SKris Buschelman 
248397b6df1SKris Buschelman     /* A-part */
249397b6df1SKris Buschelman     for (j=0; j<countA; j++){
250bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
251397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
252397b6df1SKris Buschelman       }
25316ebf90aSShri Abhyankar       val[jj++] = v1[j];
254397b6df1SKris Buschelman     }
25516ebf90aSShri Abhyankar 
25616ebf90aSShri Abhyankar     /* B-part */
25716ebf90aSShri Abhyankar     for(j=0; j < countB; j++){
258bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
259397b6df1SKris Buschelman 	row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
260397b6df1SKris Buschelman       }
26116ebf90aSShri Abhyankar       val[jj++] = v2[j];
26216ebf90aSShri Abhyankar     }
26316ebf90aSShri Abhyankar     irow++;
26416ebf90aSShri Abhyankar   }
26516ebf90aSShri Abhyankar   PetscFunctionReturn(0);
26616ebf90aSShri Abhyankar }
26716ebf90aSShri Abhyankar 
26816ebf90aSShri Abhyankar #undef __FUNCT__
26916ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpiaij"
270bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
27116ebf90aSShri Abhyankar {
27216ebf90aSShri Abhyankar   const PetscInt     *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
27316ebf90aSShri Abhyankar   PetscErrorCode     ierr;
27416ebf90aSShri Abhyankar   PetscInt           rstart,nz,i,j,jj,irow,countA,countB;
27516ebf90aSShri Abhyankar   PetscInt           *row,*col;
27616ebf90aSShri Abhyankar   const PetscScalar  *av, *bv,*v1,*v2;
27716ebf90aSShri Abhyankar   PetscScalar        *val;
27816ebf90aSShri Abhyankar   Mat_MPIAIJ         *mat =  (Mat_MPIAIJ*)A->data;
27916ebf90aSShri Abhyankar   Mat_SeqAIJ         *aa=(Mat_SeqAIJ*)(mat->A)->data;
28016ebf90aSShri Abhyankar   Mat_SeqAIJ         *bb=(Mat_SeqAIJ*)(mat->B)->data;
28116ebf90aSShri Abhyankar 
28216ebf90aSShri Abhyankar   PetscFunctionBegin;
28316ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
28416ebf90aSShri Abhyankar   garray = mat->garray;
28516ebf90aSShri Abhyankar   av=aa->a; bv=bb->a;
28616ebf90aSShri Abhyankar 
287bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX){
28816ebf90aSShri Abhyankar     nz = aa->nz + bb->nz;
28916ebf90aSShri Abhyankar     *nnz = nz;
290*185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
291*185f6596SHong Zhang     col  = row + nz;
292*185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
293*185f6596SHong Zhang 
29416ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
29516ebf90aSShri Abhyankar   } else {
29616ebf90aSShri Abhyankar     row = *r; col = *c; val = *v;
29716ebf90aSShri Abhyankar   }
29816ebf90aSShri Abhyankar 
29916ebf90aSShri Abhyankar   jj = 0; irow = rstart;
30016ebf90aSShri Abhyankar   for ( i=0; i<m; i++ ) {
30116ebf90aSShri Abhyankar     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
30216ebf90aSShri Abhyankar     countA = ai[i+1] - ai[i];
30316ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
30416ebf90aSShri Abhyankar     bjj    = bj + bi[i];
30516ebf90aSShri Abhyankar     v1     = av + ai[i];
30616ebf90aSShri Abhyankar     v2     = bv + bi[i];
30716ebf90aSShri Abhyankar 
30816ebf90aSShri Abhyankar     /* A-part */
30916ebf90aSShri Abhyankar     for (j=0; j<countA; j++){
310bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX){
31116ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
31216ebf90aSShri Abhyankar       }
31316ebf90aSShri Abhyankar       val[jj++] = v1[j];
31416ebf90aSShri Abhyankar     }
31516ebf90aSShri Abhyankar 
31616ebf90aSShri Abhyankar     /* B-part */
31716ebf90aSShri Abhyankar     for(j=0; j < countB; j++){
318bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX){
31916ebf90aSShri Abhyankar 	row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
32016ebf90aSShri Abhyankar       }
32116ebf90aSShri Abhyankar       val[jj++] = v2[j];
32216ebf90aSShri Abhyankar     }
32316ebf90aSShri Abhyankar     irow++;
32416ebf90aSShri Abhyankar   }
32516ebf90aSShri Abhyankar   PetscFunctionReturn(0);
32616ebf90aSShri Abhyankar }
32716ebf90aSShri Abhyankar 
32816ebf90aSShri Abhyankar #undef __FUNCT__
32967877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpibaij_mpiaij"
330bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
33167877ebaSShri Abhyankar {
33267877ebaSShri Abhyankar   Mat_MPIBAIJ        *mat =  (Mat_MPIBAIJ*)A->data;
33367877ebaSShri Abhyankar   Mat_SeqBAIJ        *aa=(Mat_SeqBAIJ*)(mat->A)->data;
33467877ebaSShri Abhyankar   Mat_SeqBAIJ        *bb=(Mat_SeqBAIJ*)(mat->B)->data;
33567877ebaSShri Abhyankar   const PetscInt     *ai = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj;
33667877ebaSShri Abhyankar   const PetscInt     *garray = mat->garray,mbs=mat->mbs,rstartbs=mat->rstartbs;
33767877ebaSShri Abhyankar   const PetscInt     bs = A->rmap->bs,bs2=mat->bs2;
33867877ebaSShri Abhyankar   PetscErrorCode     ierr;
33967877ebaSShri Abhyankar   PetscInt           nz,i,j,k,n,jj,irow,countA,countB,idx;
34067877ebaSShri Abhyankar   PetscInt           *row,*col;
34167877ebaSShri Abhyankar   const PetscScalar  *av=aa->a, *bv=bb->a,*v1,*v2;
34267877ebaSShri Abhyankar   PetscScalar        *val;
34367877ebaSShri Abhyankar 
34467877ebaSShri Abhyankar   PetscFunctionBegin;
34567877ebaSShri Abhyankar 
346bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
34767877ebaSShri Abhyankar     nz = bs2*(aa->nz + bb->nz);
34867877ebaSShri Abhyankar     *nnz = nz;
349*185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
350*185f6596SHong Zhang     col  = row + nz;
351*185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
352*185f6596SHong Zhang 
35367877ebaSShri Abhyankar     *r = row; *c = col; *v = val;
35467877ebaSShri Abhyankar   } else {
35567877ebaSShri Abhyankar     row = *r; col = *c; val = *v;
35667877ebaSShri Abhyankar   }
35767877ebaSShri Abhyankar 
35867877ebaSShri Abhyankar   jj = 0; irow = rstartbs;
35967877ebaSShri Abhyankar   for ( i=0; i<mbs; i++ ) {
36067877ebaSShri Abhyankar     countA = ai[i+1] - ai[i];
36167877ebaSShri Abhyankar     countB = bi[i+1] - bi[i];
36267877ebaSShri Abhyankar     ajj    = aj + ai[i];
36367877ebaSShri Abhyankar     bjj    = bj + bi[i];
36467877ebaSShri Abhyankar     v1     = av + bs2*ai[i];
36567877ebaSShri Abhyankar     v2     = bv + bs2*bi[i];
36667877ebaSShri Abhyankar 
36767877ebaSShri Abhyankar     idx = 0;
36867877ebaSShri Abhyankar     /* A-part */
36967877ebaSShri Abhyankar     for (k=0; k<countA; k++){
37067877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
37167877ebaSShri Abhyankar 	for (n=0; n<bs; n++) {
372bccb9932SShri Abhyankar 	  if (reuse == MAT_INITIAL_MATRIX){
37367877ebaSShri Abhyankar 	    row[jj] = bs*irow + n + shift;
37467877ebaSShri Abhyankar 	    col[jj] = bs*(rstartbs + ajj[k]) + j + shift;
37567877ebaSShri Abhyankar 	  }
37667877ebaSShri Abhyankar 	  val[jj++] = v1[idx++];
37767877ebaSShri Abhyankar 	}
37867877ebaSShri Abhyankar       }
37967877ebaSShri Abhyankar     }
38067877ebaSShri Abhyankar 
38167877ebaSShri Abhyankar     idx = 0;
38267877ebaSShri Abhyankar     /* B-part */
38367877ebaSShri Abhyankar     for(k=0; k<countB; k++){
38467877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
38567877ebaSShri Abhyankar 	for (n=0; n<bs; n++) {
386bccb9932SShri Abhyankar 	  if (reuse == MAT_INITIAL_MATRIX){
38767877ebaSShri Abhyankar 	    row[jj] = bs*irow + n + shift;
38867877ebaSShri Abhyankar 	    col[jj] = bs*(garray[bjj[k]]) + j + shift;
38967877ebaSShri Abhyankar 	  }
39067877ebaSShri Abhyankar 	  val[jj++] = bv[idx++];
39167877ebaSShri Abhyankar 	}
39267877ebaSShri Abhyankar       }
39367877ebaSShri Abhyankar     }
39467877ebaSShri Abhyankar     irow++;
39567877ebaSShri Abhyankar   }
39667877ebaSShri Abhyankar   PetscFunctionReturn(0);
39767877ebaSShri Abhyankar }
39867877ebaSShri Abhyankar 
39967877ebaSShri Abhyankar #undef __FUNCT__
40016ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpisbaij"
401bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
40216ebf90aSShri Abhyankar {
40316ebf90aSShri Abhyankar   const PetscInt     *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
40416ebf90aSShri Abhyankar   PetscErrorCode     ierr;
40516ebf90aSShri Abhyankar   PetscInt           rstart,nz,nza,nzb_low,i,j,jj,irow,countA,countB;
40616ebf90aSShri Abhyankar   PetscInt           *row,*col;
40716ebf90aSShri Abhyankar   const PetscScalar  *av, *bv,*v1,*v2;
40816ebf90aSShri Abhyankar   PetscScalar        *val;
40916ebf90aSShri Abhyankar   Mat_MPIAIJ         *mat =  (Mat_MPIAIJ*)A->data;
41016ebf90aSShri Abhyankar   Mat_SeqAIJ         *aa=(Mat_SeqAIJ*)(mat->A)->data;
41116ebf90aSShri Abhyankar   Mat_SeqAIJ         *bb=(Mat_SeqAIJ*)(mat->B)->data;
41216ebf90aSShri Abhyankar 
41316ebf90aSShri Abhyankar   PetscFunctionBegin;
41416ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j; adiag=aa->diag;
41516ebf90aSShri Abhyankar   bi=bb->i; bj=bb->j; garray = mat->garray;
41616ebf90aSShri Abhyankar   av=aa->a; bv=bb->a;
41716ebf90aSShri Abhyankar   rstart = A->rmap->rstart;
41816ebf90aSShri Abhyankar 
419bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
42016ebf90aSShri Abhyankar     nza = 0;nzb_low = 0;
42116ebf90aSShri Abhyankar     for(i=0; i<m; i++){
42216ebf90aSShri Abhyankar       nza     = nza + (ai[i+1] - adiag[i]);
42316ebf90aSShri Abhyankar       countB  = bi[i+1] - bi[i];
42416ebf90aSShri Abhyankar       bjj     = bj + bi[i];
42516ebf90aSShri Abhyankar 
42616ebf90aSShri Abhyankar       j = 0;
42716ebf90aSShri Abhyankar       while(garray[bjj[j]] < rstart) {
42816ebf90aSShri Abhyankar 	if(j == countB) break;
42916ebf90aSShri Abhyankar 	j++;nzb_low++;
43016ebf90aSShri Abhyankar       }
43116ebf90aSShri Abhyankar     }
43216ebf90aSShri Abhyankar     /* Total nz = nz for the upper triangular A part + nz for the 2nd B part */
43316ebf90aSShri Abhyankar     nz = nza + (bb->nz - nzb_low);
43416ebf90aSShri Abhyankar     *nnz = nz;
435*185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
436*185f6596SHong Zhang     col  = row + nz;
437*185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
438*185f6596SHong Zhang 
43916ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
44016ebf90aSShri Abhyankar   } else {
44116ebf90aSShri Abhyankar     row = *r; col = *c; val = *v;
44216ebf90aSShri Abhyankar   }
44316ebf90aSShri Abhyankar 
44416ebf90aSShri Abhyankar   jj = 0; irow = rstart;
44516ebf90aSShri Abhyankar   for ( i=0; i<m; i++ ) {
44616ebf90aSShri Abhyankar     ajj    = aj + adiag[i];                 /* ptr to the beginning of the diagonal of this row */
44716ebf90aSShri Abhyankar     v1     = av + adiag[i];
44816ebf90aSShri Abhyankar     countA = ai[i+1] - adiag[i];
44916ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
45016ebf90aSShri Abhyankar     bjj    = bj + bi[i];
45116ebf90aSShri Abhyankar     v2     = bv + bi[i];
45216ebf90aSShri Abhyankar 
45316ebf90aSShri Abhyankar      /* A-part */
45416ebf90aSShri Abhyankar     for (j=0; j<countA; j++){
455bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
45616ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
45716ebf90aSShri Abhyankar       }
45816ebf90aSShri Abhyankar       val[jj++] = v1[j];
45916ebf90aSShri Abhyankar     }
46016ebf90aSShri Abhyankar 
46116ebf90aSShri Abhyankar     /* B-part */
46216ebf90aSShri Abhyankar     for(j=0; j < countB; j++){
46316ebf90aSShri Abhyankar       if (garray[bjj[j]] > rstart) {
464bccb9932SShri Abhyankar 	if (reuse == MAT_INITIAL_MATRIX) {
46516ebf90aSShri Abhyankar 	  row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
46616ebf90aSShri Abhyankar 	}
46716ebf90aSShri Abhyankar 	val[jj++] = v2[j];
46816ebf90aSShri Abhyankar       }
469397b6df1SKris Buschelman     }
470397b6df1SKris Buschelman     irow++;
471397b6df1SKris Buschelman   }
472397b6df1SKris Buschelman   PetscFunctionReturn(0);
473397b6df1SKris Buschelman }
474397b6df1SKris Buschelman 
475397b6df1SKris Buschelman #undef __FUNCT__
4763924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS"
477dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A)
478dfbe8321SBarry Smith {
479f0c56d0fSKris Buschelman   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
480dfbe8321SBarry Smith   PetscErrorCode ierr;
481c1490034SHong Zhang   PetscMPIInt    size=lu->size;
482b24902e0SBarry Smith 
483397b6df1SKris Buschelman   PetscFunctionBegin;
484397b6df1SKris Buschelman   if (lu->CleanUpMUMPS) {
485397b6df1SKris Buschelman     /* Terminate instance, deallocate memories */
486329ec9b3SHong Zhang     if (size > 1){
48768653410SBarry Smith       ierr = PetscFree2(lu->id.sol_loc,lu->id.isol_loc);CHKERRQ(ierr);
488329ec9b3SHong Zhang       ierr = VecScatterDestroy(lu->scat_rhs);CHKERRQ(ierr);
489329ec9b3SHong Zhang       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);}
492450b117fSShri Abhyankar     }
493*185f6596SHong 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);
50497969023SHong Zhang   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   /*-------------*/
543397b6df1SKris Buschelman   lu->id.job = 3;
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       }
557329ec9b3SHong Zhang       ierr = ISCreateGeneral(PETSC_COMM_SELF,lu->id.lsol_loc,lu->id.isol_loc,&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 
569ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX)
570a58c3f20SHong Zhang /*
571a58c3f20SHong Zhang   input:
572a58c3f20SHong Zhang    F:        numeric factor
573a58c3f20SHong Zhang   output:
574a58c3f20SHong Zhang    nneg:     total number of negative pivots
575a58c3f20SHong Zhang    nzero:    0
576a58c3f20SHong Zhang    npos:     (global dimension of F) - nneg
577a58c3f20SHong Zhang */
578a58c3f20SHong Zhang 
579a58c3f20SHong Zhang #undef __FUNCT__
580a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
581dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
582a58c3f20SHong Zhang {
583a58c3f20SHong Zhang   Mat_MUMPS      *lu =(Mat_MUMPS*)F->spptr;
584dfbe8321SBarry Smith   PetscErrorCode ierr;
585c1490034SHong Zhang   PetscMPIInt    size;
586a58c3f20SHong Zhang 
587a58c3f20SHong Zhang   PetscFunctionBegin;
5887adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)F)->comm,&size);CHKERRQ(ierr);
589bcb30aebSHong 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 */
59065e19b50SBarry 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));
591a58c3f20SHong Zhang   if (nneg){
592a58c3f20SHong Zhang     if (!lu->myid){
593a58c3f20SHong Zhang       *nneg = lu->id.INFOG(12);
594a58c3f20SHong Zhang     }
595bcb30aebSHong Zhang     ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr);
596a58c3f20SHong Zhang   }
597a58c3f20SHong Zhang   if (nzero) *nzero = 0;
598d0f46423SBarry Smith   if (npos)  *npos  = F->rmap->N - (*nneg);
599a58c3f20SHong Zhang   PetscFunctionReturn(0);
600a58c3f20SHong Zhang }
601ace3df97SHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */
602a58c3f20SHong Zhang 
603397b6df1SKris Buschelman #undef __FUNCT__
604f6c57405SHong Zhang #define __FUNCT__ "MatFactorNumeric_MUMPS"
6050481f469SBarry Smith PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
606af281ebdSHong Zhang {
607dcd589f8SShri Abhyankar   Mat_MUMPS       *lu =(Mat_MUMPS*)(F)->spptr;
6086849ba73SBarry Smith   PetscErrorCode  ierr;
609bccb9932SShri Abhyankar   MatReuse        reuse;
610e09efc27SHong Zhang   Mat             F_diag;
61116ebf90aSShri Abhyankar   PetscTruth      isMPIAIJ;
612397b6df1SKris Buschelman 
613397b6df1SKris Buschelman   PetscFunctionBegin;
614bccb9932SShri Abhyankar   reuse = MAT_REUSE_MATRIX;
615bccb9932SShri Abhyankar   ierr = (*lu->ConvertToTriples)(A, 1, reuse, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
616397b6df1SKris Buschelman 
617397b6df1SKris Buschelman   /* numerical factorization phase */
618329ec9b3SHong Zhang   /*-------------------------------*/
619329ec9b3SHong Zhang   lu->id.job = 2;
620958c9bccSBarry Smith   if(!lu->id.ICNTL(18)) {
621a7aca84bSHong Zhang     if (!lu->myid) {
622397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
623397b6df1SKris Buschelman       lu->id.a = (mumps_double_complex*)lu->val;
624397b6df1SKris Buschelman #else
625397b6df1SKris Buschelman       lu->id.a = lu->val;
626397b6df1SKris Buschelman #endif
627397b6df1SKris Buschelman     }
628397b6df1SKris Buschelman   } else {
629397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
630397b6df1SKris Buschelman     lu->id.a_loc = (mumps_double_complex*)lu->val;
631397b6df1SKris Buschelman #else
632397b6df1SKris Buschelman     lu->id.a_loc = lu->val;
633397b6df1SKris Buschelman #endif
634397b6df1SKris Buschelman   }
635397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
636397b6df1SKris Buschelman   zmumps_c(&lu->id);
637397b6df1SKris Buschelman #else
638397b6df1SKris Buschelman   dmumps_c(&lu->id);
639397b6df1SKris Buschelman #endif
640397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
64165e19b50SBarry 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));
64265e19b50SBarry 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));
643397b6df1SKris Buschelman   }
64465e19b50SBarry 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));
645397b6df1SKris Buschelman 
6468ada1bb4SHong Zhang   if (lu->size > 1){
64716ebf90aSShri Abhyankar     ierr = PetscTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr);
648a214ac2aSShri Abhyankar     if(isMPIAIJ) {
649dcd589f8SShri Abhyankar       F_diag = ((Mat_MPIAIJ *)(F)->data)->A;
650e09efc27SHong Zhang     } else {
651dcd589f8SShri Abhyankar       F_diag = ((Mat_MPISBAIJ *)(F)->data)->A;
652e09efc27SHong Zhang     }
653e09efc27SHong Zhang     F_diag->assembled = PETSC_TRUE;
654329ec9b3SHong Zhang     if (lu->nSolve){
655329ec9b3SHong Zhang       ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);
6560e83c824SBarry Smith       ierr = PetscFree2(lu->id.sol_loc,lu->id.isol_loc);CHKERRQ(ierr);
657329ec9b3SHong Zhang       ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);
658329ec9b3SHong Zhang     }
6598ada1bb4SHong Zhang   }
660dcd589f8SShri Abhyankar   (F)->assembled   = PETSC_TRUE;
661397b6df1SKris Buschelman   lu->matstruc     = SAME_NONZERO_PATTERN;
662ace87b0dSHong Zhang   lu->CleanUpMUMPS = PETSC_TRUE;
663329ec9b3SHong Zhang   lu->nSolve       = 0;
66467877ebaSShri Abhyankar 
66567877ebaSShri Abhyankar   if (lu->size > 1){
66667877ebaSShri Abhyankar     /* distributed solution */
66767877ebaSShri Abhyankar     lu->id.ICNTL(21) = 1;
66867877ebaSShri Abhyankar     if (!lu->nSolve){
66967877ebaSShri Abhyankar       /* Create x_seq=sol_loc for repeated use */
67067877ebaSShri Abhyankar       PetscInt    lsol_loc;
67167877ebaSShri Abhyankar       PetscScalar *sol_loc;
67267877ebaSShri Abhyankar       lsol_loc = lu->id.INFO(23); /* length of sol_loc */
67367877ebaSShri Abhyankar       ierr = PetscMalloc2(lsol_loc,PetscScalar,&sol_loc,lsol_loc,PetscInt,&lu->id.isol_loc);CHKERRQ(ierr);
67467877ebaSShri Abhyankar       lu->id.lsol_loc = lsol_loc;
67567877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
67667877ebaSShri Abhyankar       lu->id.sol_loc  = (mumps_double_complex*)sol_loc;
67767877ebaSShri Abhyankar #else
67867877ebaSShri Abhyankar       lu->id.sol_loc  = sol_loc;
67967877ebaSShri Abhyankar #endif
68067877ebaSShri Abhyankar       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,lsol_loc,sol_loc,&lu->x_seq);CHKERRQ(ierr);
68167877ebaSShri Abhyankar     }
68267877ebaSShri Abhyankar   }
683397b6df1SKris Buschelman   PetscFunctionReturn(0);
684397b6df1SKris Buschelman }
685397b6df1SKris Buschelman 
686dcd589f8SShri Abhyankar #undef __FUNCT__
687dcd589f8SShri Abhyankar #define __FUNCT__ "PetscSetMUMPSOptions"
688dcd589f8SShri Abhyankar PetscErrorCode PetscSetMUMPSOptions(Mat F, Mat A)
689dcd589f8SShri Abhyankar {
690dcd589f8SShri Abhyankar   Mat_MUMPS        *lu = (Mat_MUMPS*)F->spptr;
691dcd589f8SShri Abhyankar   PetscErrorCode   ierr;
692dcd589f8SShri Abhyankar   PetscInt         icntl;
693dcd589f8SShri Abhyankar   PetscTruth       flg;
694dcd589f8SShri Abhyankar 
695dcd589f8SShri Abhyankar   PetscFunctionBegin;
696dcd589f8SShri Abhyankar   ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
697cb828f0fSHong Zhang   ierr = PetscOptionsTruth("-mat_mumps_view","View MUMPS parameters","None",lu->mumpsview,&lu->mumpsview,PETSC_NULL);CHKERRQ(ierr);
698dcd589f8SShri Abhyankar   if (lu->size == 1){
699dcd589f8SShri Abhyankar     lu->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
700dcd589f8SShri Abhyankar   } else {
701dcd589f8SShri Abhyankar     lu->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
702dcd589f8SShri Abhyankar   }
703dcd589f8SShri Abhyankar 
704dcd589f8SShri Abhyankar   icntl=-1;
705dcd589f8SShri Abhyankar   lu->id.ICNTL(4) = 0;  /* level of printing; overwrite mumps default ICNTL(4)=2 */
706dcd589f8SShri Abhyankar   ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
707dcd589f8SShri Abhyankar   if ((flg && icntl > 0) || PetscLogPrintInfo) {
708dcd589f8SShri Abhyankar     lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */
709dcd589f8SShri Abhyankar   } else { /* no output */
710dcd589f8SShri Abhyankar     lu->id.ICNTL(1) = 0;  /* error message, default= 6 */
711dcd589f8SShri Abhyankar     lu->id.ICNTL(2) = 0;  /* output stream for diagnostic printing, statistics, and warning. default=0 */
712dcd589f8SShri Abhyankar     lu->id.ICNTL(3) = 0; /* output stream for global information, default=6 */
713dcd589f8SShri Abhyankar   }
714dcd589f8SShri 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);
715dcd589f8SShri Abhyankar   icntl=-1;
716dcd589f8SShri Abhyankar   ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
717dcd589f8SShri Abhyankar   if (flg) {
718dcd589f8SShri Abhyankar     if (icntl== 1){
719e32f2f54SBarry 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");
720dcd589f8SShri Abhyankar     } else {
721dcd589f8SShri Abhyankar       lu->id.ICNTL(7) = icntl;
722dcd589f8SShri Abhyankar     }
723dcd589f8SShri Abhyankar   }
724dcd589f8SShri 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);
725dcd589f8SShri Abhyankar   ierr = PetscOptionsInt("-mat_mumps_icntl_9","ICNTL(9): A or A^T x=b to be solved. 1: A; otherwise: A^T","None",lu->id.ICNTL(9),&lu->id.ICNTL(9),PETSC_NULL);CHKERRQ(ierr);
726dcd589f8SShri 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);
727dcd589f8SShri 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);
728dcd589f8SShri 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);
729dcd589f8SShri 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);
730dcd589f8SShri 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);
731dcd589f8SShri Abhyankar   ierr = PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): Schur complement","None",lu->id.ICNTL(19),&lu->id.ICNTL(19),PETSC_NULL);CHKERRQ(ierr);
732dcd589f8SShri Abhyankar 
733dcd589f8SShri 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);
734dcd589f8SShri 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);
735dcd589f8SShri 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);
736dcd589f8SShri 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);
737dcd589f8SShri 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);
738dcd589f8SShri Abhyankar   ierr = PetscOptionsInt("-mat_mumps_icntl_27","ICNTL(27): experimental parameter","None",lu->id.ICNTL(27),&lu->id.ICNTL(27),PETSC_NULL);CHKERRQ(ierr);
739dcd589f8SShri Abhyankar 
740dcd589f8SShri 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);
741dcd589f8SShri 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);
742dcd589f8SShri 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);
743dcd589f8SShri 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);
744dcd589f8SShri 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);
745dcd589f8SShri Abhyankar   PetscOptionsEnd();
746dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
747dcd589f8SShri Abhyankar }
748dcd589f8SShri Abhyankar 
749dcd589f8SShri Abhyankar #undef __FUNCT__
750dcd589f8SShri Abhyankar #define __FUNCT__ "PetscInitializeMUMPS"
751dcd589f8SShri Abhyankar PetscErrorCode PetscInitializeMUMPS(Mat F)
752dcd589f8SShri Abhyankar {
753dcd589f8SShri Abhyankar   Mat_MUMPS       *lu = (Mat_MUMPS*)F->spptr;
754dcd589f8SShri Abhyankar   PetscErrorCode  ierr;
755dcd589f8SShri Abhyankar   PetscInt        icntl;
756dcd589f8SShri Abhyankar   PetscTruth      flg;
757dcd589f8SShri Abhyankar 
758dcd589f8SShri Abhyankar   PetscFunctionBegin;
759dcd589f8SShri Abhyankar   lu->id.job = JOB_INIT;
760dcd589f8SShri Abhyankar   lu->id.par=1;  /* host participates factorizaton and solve */
761dcd589f8SShri Abhyankar   lu->id.sym=lu->sym;
762dcd589f8SShri Abhyankar   if (lu->sym == 2){
763dcd589f8SShri Abhyankar     ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr);
764dcd589f8SShri Abhyankar     if (flg && icntl == 1) lu->id.sym=icntl;  /* matrix is spd */
765dcd589f8SShri Abhyankar   }
766dcd589f8SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
767dcd589f8SShri Abhyankar   zmumps_c(&lu->id);
768dcd589f8SShri Abhyankar #else
769dcd589f8SShri Abhyankar   dmumps_c(&lu->id);
770dcd589f8SShri Abhyankar #endif
771dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
772dcd589f8SShri Abhyankar }
773dcd589f8SShri Abhyankar 
774397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */
775397b6df1SKris Buschelman #undef __FUNCT__
776f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
7770481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
778b24902e0SBarry Smith {
779719d5645SBarry Smith   Mat_MUMPS          *lu = (Mat_MUMPS*)F->spptr;
780dcd589f8SShri Abhyankar   PetscErrorCode     ierr;
781bccb9932SShri Abhyankar   MatReuse           reuse;
78267877ebaSShri Abhyankar   Vec                b;
78367877ebaSShri Abhyankar   IS                 is_iden;
78467877ebaSShri Abhyankar   const PetscInt     M = A->rmap->N;
785397b6df1SKris Buschelman 
786397b6df1SKris Buschelman   PetscFunctionBegin;
787b24902e0SBarry Smith   lu->sym                  = 0;
788b24902e0SBarry Smith   lu->matstruc             = DIFFERENT_NONZERO_PATTERN;
789dcd589f8SShri Abhyankar 
790dcd589f8SShri Abhyankar   ierr = MPI_Comm_rank(((PetscObject)A)->comm, &lu->myid);
791dcd589f8SShri Abhyankar   ierr = MPI_Comm_size(((PetscObject)A)->comm,&lu->size);CHKERRQ(ierr);
792dcd589f8SShri Abhyankar   ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->comm_mumps));CHKERRQ(ierr);
793dcd589f8SShri Abhyankar   lu->id.comm_fortran = MPI_Comm_c2f(lu->comm_mumps);
794dcd589f8SShri Abhyankar 
795dcd589f8SShri Abhyankar   /* Initialize a MUMPS instance */
796dcd589f8SShri Abhyankar   ierr = PetscInitializeMUMPS(F);CHKERRQ(ierr);
797dcd589f8SShri Abhyankar   /* Set MUMPS options */
798dcd589f8SShri Abhyankar   ierr = PetscSetMUMPSOptions(F,A);CHKERRQ(ierr);
799dcd589f8SShri Abhyankar 
800bccb9932SShri Abhyankar   reuse = MAT_INITIAL_MATRIX;
801bccb9932SShri Abhyankar   ierr = (*lu->ConvertToTriples)(A, 1, reuse, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
802dcd589f8SShri Abhyankar 
80367877ebaSShri Abhyankar   /* analysis phase */
80467877ebaSShri Abhyankar   /*----------------*/
80567877ebaSShri Abhyankar   lu->id.job = 1;
80667877ebaSShri Abhyankar   lu->id.n = M;
80767877ebaSShri Abhyankar   switch (lu->id.ICNTL(18)){
80867877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
80967877ebaSShri Abhyankar     if (!lu->myid) {
81067877ebaSShri Abhyankar       lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
81167877ebaSShri Abhyankar       if (lu->id.ICNTL(6)>1){
81267877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
81367877ebaSShri Abhyankar         lu->id.a = (mumps_double_complex*)lu->val;
81467877ebaSShri Abhyankar #else
81567877ebaSShri Abhyankar         lu->id.a = lu->val;
81667877ebaSShri Abhyankar #endif
81767877ebaSShri Abhyankar       }
81867877ebaSShri Abhyankar     }
81967877ebaSShri Abhyankar     break;
82067877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
82167877ebaSShri Abhyankar     lu->id.nz_loc = lu->nz;
82267877ebaSShri Abhyankar     lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
82367877ebaSShri Abhyankar     if (lu->id.ICNTL(6)>1) {
82467877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
82567877ebaSShri Abhyankar       lu->id.a_loc = (mumps_double_complex*)lu->val;
82667877ebaSShri Abhyankar #else
82767877ebaSShri Abhyankar       lu->id.a_loc = lu->val;
82867877ebaSShri Abhyankar #endif
82967877ebaSShri Abhyankar     }
83067877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
83167877ebaSShri Abhyankar     if (!lu->myid){
83267877ebaSShri Abhyankar       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr);
83367877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
83467877ebaSShri Abhyankar     } else {
83567877ebaSShri Abhyankar       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr);
83667877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
83767877ebaSShri Abhyankar     }
83867877ebaSShri Abhyankar     ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr);
83967877ebaSShri Abhyankar     ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
84067877ebaSShri Abhyankar     ierr = VecSetFromOptions(b);CHKERRQ(ierr);
84167877ebaSShri Abhyankar 
84267877ebaSShri Abhyankar     ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr);
84367877ebaSShri Abhyankar     ierr = ISDestroy(is_iden);CHKERRQ(ierr);
84467877ebaSShri Abhyankar     ierr = VecDestroy(b);CHKERRQ(ierr);
84567877ebaSShri Abhyankar     break;
84667877ebaSShri Abhyankar     }
84767877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
84867877ebaSShri Abhyankar   zmumps_c(&lu->id);
84967877ebaSShri Abhyankar #else
85067877ebaSShri Abhyankar   dmumps_c(&lu->id);
85167877ebaSShri Abhyankar #endif
85267877ebaSShri 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));
85367877ebaSShri Abhyankar 
854719d5645SBarry Smith   F->ops->lufactornumeric  = MatFactorNumeric_MUMPS;
855dcd589f8SShri Abhyankar   F->ops->solve            = MatSolve_MUMPS;
856b24902e0SBarry Smith   PetscFunctionReturn(0);
857b24902e0SBarry Smith }
858b24902e0SBarry Smith 
859450b117fSShri Abhyankar /* Note the Petsc r and c permutations are ignored */
860450b117fSShri Abhyankar #undef __FUNCT__
861450b117fSShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_BAIJMUMPS"
862450b117fSShri Abhyankar PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
863450b117fSShri Abhyankar {
864dcd589f8SShri Abhyankar 
865450b117fSShri Abhyankar   Mat_MUMPS       *lu = (Mat_MUMPS*)F->spptr;
866dcd589f8SShri Abhyankar   PetscErrorCode  ierr;
867bccb9932SShri Abhyankar   MatReuse        reuse;
86867877ebaSShri Abhyankar   Vec             b;
86967877ebaSShri Abhyankar   IS              is_iden;
87067877ebaSShri Abhyankar   const PetscInt  M = A->rmap->N;
871450b117fSShri Abhyankar 
872450b117fSShri Abhyankar   PetscFunctionBegin;
873450b117fSShri Abhyankar   lu->sym                  = 0;
874450b117fSShri Abhyankar   lu->matstruc             = DIFFERENT_NONZERO_PATTERN;
875dcd589f8SShri Abhyankar   ierr = MPI_Comm_rank(((PetscObject)A)->comm, &lu->myid);
876dcd589f8SShri Abhyankar   ierr = MPI_Comm_size(((PetscObject)A)->comm,&lu->size);CHKERRQ(ierr);
877dcd589f8SShri Abhyankar   ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->comm_mumps));CHKERRQ(ierr);
878dcd589f8SShri Abhyankar   lu->id.comm_fortran = MPI_Comm_c2f(lu->comm_mumps);
879dcd589f8SShri Abhyankar 
880dcd589f8SShri Abhyankar   /* Initialize a MUMPS instance */
881dcd589f8SShri Abhyankar   ierr = PetscInitializeMUMPS(F);CHKERRQ(ierr);
882dcd589f8SShri Abhyankar   /* Set MUMPS options */
883dcd589f8SShri Abhyankar   ierr = PetscSetMUMPSOptions(F,A);CHKERRQ(ierr);
884dcd589f8SShri Abhyankar 
885bccb9932SShri Abhyankar   reuse = MAT_INITIAL_MATRIX;
886bccb9932SShri Abhyankar   ierr = (*lu->ConvertToTriples)(A, 1, reuse, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
88767877ebaSShri Abhyankar 
88867877ebaSShri Abhyankar   /* analysis phase */
88967877ebaSShri Abhyankar   /*----------------*/
89067877ebaSShri Abhyankar   lu->id.job = 1;
89167877ebaSShri Abhyankar   lu->id.n = M;
89267877ebaSShri Abhyankar   switch (lu->id.ICNTL(18)){
89367877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
89467877ebaSShri Abhyankar     if (!lu->myid) {
89567877ebaSShri Abhyankar       lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
89667877ebaSShri Abhyankar       if (lu->id.ICNTL(6)>1){
89767877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
89867877ebaSShri Abhyankar         lu->id.a = (mumps_double_complex*)lu->val;
89967877ebaSShri Abhyankar #else
90067877ebaSShri Abhyankar         lu->id.a = lu->val;
90167877ebaSShri Abhyankar #endif
90267877ebaSShri Abhyankar       }
90367877ebaSShri Abhyankar     }
90467877ebaSShri Abhyankar     break;
90567877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
90667877ebaSShri Abhyankar     lu->id.nz_loc = lu->nz;
90767877ebaSShri Abhyankar     lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
90867877ebaSShri Abhyankar     if (lu->id.ICNTL(6)>1) {
90967877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
91067877ebaSShri Abhyankar       lu->id.a_loc = (mumps_double_complex*)lu->val;
91167877ebaSShri Abhyankar #else
91267877ebaSShri Abhyankar       lu->id.a_loc = lu->val;
91367877ebaSShri Abhyankar #endif
91467877ebaSShri Abhyankar     }
91567877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
91667877ebaSShri Abhyankar     if (!lu->myid){
91767877ebaSShri Abhyankar       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr);
91867877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
91967877ebaSShri Abhyankar     } else {
92067877ebaSShri Abhyankar       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr);
92167877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
92267877ebaSShri Abhyankar     }
92367877ebaSShri Abhyankar     ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr);
92467877ebaSShri Abhyankar     ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
92567877ebaSShri Abhyankar     ierr = VecSetFromOptions(b);CHKERRQ(ierr);
92667877ebaSShri Abhyankar 
92767877ebaSShri Abhyankar     ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr);
92867877ebaSShri Abhyankar     ierr = ISDestroy(is_iden);CHKERRQ(ierr);
92967877ebaSShri Abhyankar     ierr = VecDestroy(b);CHKERRQ(ierr);
93067877ebaSShri Abhyankar     break;
93167877ebaSShri Abhyankar     }
93267877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
93367877ebaSShri Abhyankar   zmumps_c(&lu->id);
93467877ebaSShri Abhyankar #else
93567877ebaSShri Abhyankar   dmumps_c(&lu->id);
93667877ebaSShri Abhyankar #endif
93767877ebaSShri 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));
93867877ebaSShri Abhyankar 
939450b117fSShri Abhyankar   F->ops->lufactornumeric  = MatFactorNumeric_MUMPS;
940dcd589f8SShri Abhyankar   F->ops->solve            = MatSolve_MUMPS;
941450b117fSShri Abhyankar   PetscFunctionReturn(0);
942450b117fSShri Abhyankar }
943b24902e0SBarry Smith 
944397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */
945397b6df1SKris Buschelman #undef __FUNCT__
94667877ebaSShri Abhyankar #define __FUNCT__ "MatCholeskyFactorSymbolic_MUMPS"
94767877ebaSShri Abhyankar PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
948b24902e0SBarry Smith {
9492792810eSHong Zhang   Mat_MUMPS          *lu = (Mat_MUMPS*)F->spptr;
950dcd589f8SShri Abhyankar   PetscErrorCode     ierr;
951bccb9932SShri Abhyankar   MatReuse           reuse;
95267877ebaSShri Abhyankar   Vec                b;
95367877ebaSShri Abhyankar   IS                 is_iden;
95467877ebaSShri Abhyankar   const PetscInt     M = A->rmap->N;
955397b6df1SKris Buschelman 
956397b6df1SKris Buschelman   PetscFunctionBegin;
957b24902e0SBarry Smith   lu->sym                          = 2;
958b24902e0SBarry Smith   lu->matstruc                     = DIFFERENT_NONZERO_PATTERN;
959dcd589f8SShri Abhyankar   ierr = MPI_Comm_rank(((PetscObject)A)->comm, &lu->myid);
960dcd589f8SShri Abhyankar   ierr = MPI_Comm_size(((PetscObject)A)->comm,&lu->size);CHKERRQ(ierr);
961dcd589f8SShri Abhyankar   ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->comm_mumps));CHKERRQ(ierr);
962dcd589f8SShri Abhyankar   lu->id.comm_fortran = MPI_Comm_c2f(lu->comm_mumps);
963dcd589f8SShri Abhyankar 
964dcd589f8SShri Abhyankar   /* Initialize a MUMPS instance */
965dcd589f8SShri Abhyankar   ierr = PetscInitializeMUMPS(F);CHKERRQ(ierr);
966dcd589f8SShri Abhyankar   /* Set MUMPS options */
967dcd589f8SShri Abhyankar   ierr = PetscSetMUMPSOptions(F,A);CHKERRQ(ierr);
968dcd589f8SShri Abhyankar 
969bccb9932SShri Abhyankar   reuse = MAT_INITIAL_MATRIX;
970bccb9932SShri Abhyankar   ierr = (*lu->ConvertToTriples)(A, 1 , reuse, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
971dcd589f8SShri Abhyankar 
97267877ebaSShri Abhyankar   /* analysis phase */
97367877ebaSShri Abhyankar   /*----------------*/
97467877ebaSShri Abhyankar   lu->id.job = 1;
97567877ebaSShri Abhyankar   lu->id.n = M;
97667877ebaSShri Abhyankar   switch (lu->id.ICNTL(18)){
97767877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
97867877ebaSShri Abhyankar     if (!lu->myid) {
97967877ebaSShri Abhyankar       lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
98067877ebaSShri Abhyankar       if (lu->id.ICNTL(6)>1){
98167877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
98267877ebaSShri Abhyankar         lu->id.a = (mumps_double_complex*)lu->val;
98367877ebaSShri Abhyankar #else
98467877ebaSShri Abhyankar         lu->id.a = lu->val;
98567877ebaSShri Abhyankar #endif
98667877ebaSShri Abhyankar       }
98767877ebaSShri Abhyankar     }
98867877ebaSShri Abhyankar     break;
98967877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
99067877ebaSShri Abhyankar     lu->id.nz_loc = lu->nz;
99167877ebaSShri Abhyankar     lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
99267877ebaSShri Abhyankar     if (lu->id.ICNTL(6)>1) {
99367877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
99467877ebaSShri Abhyankar       lu->id.a_loc = (mumps_double_complex*)lu->val;
99567877ebaSShri Abhyankar #else
99667877ebaSShri Abhyankar       lu->id.a_loc = lu->val;
99767877ebaSShri Abhyankar #endif
99867877ebaSShri Abhyankar     }
99967877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
100067877ebaSShri Abhyankar     if (!lu->myid){
100167877ebaSShri Abhyankar       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr);
100267877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
100367877ebaSShri Abhyankar     } else {
100467877ebaSShri Abhyankar       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr);
100567877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
100667877ebaSShri Abhyankar     }
100767877ebaSShri Abhyankar     ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr);
100867877ebaSShri Abhyankar     ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
100967877ebaSShri Abhyankar     ierr = VecSetFromOptions(b);CHKERRQ(ierr);
101067877ebaSShri Abhyankar 
101167877ebaSShri Abhyankar     ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr);
101267877ebaSShri Abhyankar     ierr = ISDestroy(is_iden);CHKERRQ(ierr);
101367877ebaSShri Abhyankar     ierr = VecDestroy(b);CHKERRQ(ierr);
101467877ebaSShri Abhyankar     break;
101567877ebaSShri Abhyankar     }
101667877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
101767877ebaSShri Abhyankar   zmumps_c(&lu->id);
101867877ebaSShri Abhyankar #else
101967877ebaSShri Abhyankar   dmumps_c(&lu->id);
102067877ebaSShri Abhyankar #endif
102167877ebaSShri 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));
102267877ebaSShri Abhyankar 
10232792810eSHong Zhang   F->ops->choleskyfactornumeric =  MatFactorNumeric_MUMPS;
1024dcd589f8SShri Abhyankar   F->ops->solve                 =  MatSolve_MUMPS;
1025db4efbfdSBarry Smith #if !defined(PETSC_USE_COMPLEX)
1026dcd589f8SShri Abhyankar   (F)->ops->getinertia          =  MatGetInertia_SBAIJMUMPS;
1027db4efbfdSBarry Smith #endif
1028b24902e0SBarry Smith   PetscFunctionReturn(0);
1029b24902e0SBarry Smith }
1030b24902e0SBarry Smith 
1031397b6df1SKris Buschelman #undef __FUNCT__
1032f6c57405SHong Zhang #define __FUNCT__ "MatFactorInfo_MUMPS"
103374ed9c26SBarry Smith PetscErrorCode MatFactorInfo_MUMPS(Mat A,PetscViewer viewer)
103474ed9c26SBarry Smith {
1035f6c57405SHong Zhang   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
1036f6c57405SHong Zhang   PetscErrorCode ierr;
1037f6c57405SHong Zhang 
1038f6c57405SHong Zhang   PetscFunctionBegin;
1039f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):        %d \n",lu->id.ICNTL(1));CHKERRQ(ierr);
1040f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg):%d \n",lu->id.ICNTL(2));CHKERRQ(ierr);
1041f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):  %d \n",lu->id.ICNTL(3));CHKERRQ(ierr);
1042f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):       %d \n",lu->id.ICNTL(4));CHKERRQ(ierr);
1043f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):        %d \n",lu->id.ICNTL(5));CHKERRQ(ierr);
1044f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):       %d \n",lu->id.ICNTL(6));CHKERRQ(ierr);
1045f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (matrix ordering):         %d \n",lu->id.ICNTL(7));CHKERRQ(ierr);
1046f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scalling strategy):       %d \n",lu->id.ICNTL(8));CHKERRQ(ierr);
1047f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(9) (A/A^T x=b is solved):     %d \n",lu->id.ICNTL(9));CHKERRQ(ierr);
1048f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr);
1049f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):         %d \n",lu->id.ICNTL(11));CHKERRQ(ierr);
1050f6c57405SHong Zhang   if (!lu->myid && lu->id.ICNTL(11)>0) {
1051f6c57405SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(4) (inf norm of input mat):        %g\n",lu->id.RINFOG(4));CHKERRQ(ierr);
1052f6c57405SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(5) (inf norm of solution):         %g\n",lu->id.RINFOG(5));CHKERRQ(ierr);
1053f6c57405SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(6) (inf norm of residual):         %g\n",lu->id.RINFOG(6));CHKERRQ(ierr);
1054f6c57405SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",lu->id.RINFOG(7),lu->id.RINFOG(8));CHKERRQ(ierr);
1055f6c57405SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(9) (error estimate):               %g \n",lu->id.RINFOG(9));CHKERRQ(ierr);
1056f6c57405SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",lu->id.RINFOG(10),lu->id.RINFOG(11));CHKERRQ(ierr);
1057f6c57405SHong Zhang 
1058f6c57405SHong Zhang   }
1059f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",lu->id.ICNTL(12));CHKERRQ(ierr);
1060f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",lu->id.ICNTL(13));CHKERRQ(ierr);
1061f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr);
1062f6c57405SHong Zhang   /* ICNTL(15-17) not used */
1063f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",lu->id.ICNTL(18));CHKERRQ(ierr);
1064f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Shur complement info):                       %d \n",lu->id.ICNTL(19));CHKERRQ(ierr);
1065f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",lu->id.ICNTL(20));CHKERRQ(ierr);
1066f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d \n",lu->id.ICNTL(21));CHKERRQ(ierr);
1067c0165424SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",lu->id.ICNTL(22));CHKERRQ(ierr);
1068c0165424SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",lu->id.ICNTL(23));CHKERRQ(ierr);
1069c0165424SHong Zhang 
1070c0165424SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",lu->id.ICNTL(24));CHKERRQ(ierr);
1071c0165424SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",lu->id.ICNTL(25));CHKERRQ(ierr);
1072c0165424SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",lu->id.ICNTL(26));CHKERRQ(ierr);
1073c0165424SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",lu->id.ICNTL(27));CHKERRQ(ierr);
1074f6c57405SHong Zhang 
1075f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",lu->id.CNTL(1));CHKERRQ(ierr);
1076f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr);
1077f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",lu->id.CNTL(3));CHKERRQ(ierr);
1078f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",lu->id.CNTL(4));CHKERRQ(ierr);
1079c0165424SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",lu->id.CNTL(5));CHKERRQ(ierr);
1080f6c57405SHong Zhang 
1081f6c57405SHong Zhang   /* infomation local to each processor */
1082f6c57405SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "              RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);}
10837adad957SLisandro Dalcin   ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm,"              [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr);
10847adad957SLisandro Dalcin   ierr = PetscSynchronizedFlush(((PetscObject)A)->comm);
1085f6c57405SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "              RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);}
10867adad957SLisandro Dalcin   ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm,"              [%d]  %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr);
10877adad957SLisandro Dalcin   ierr = PetscSynchronizedFlush(((PetscObject)A)->comm);
1088f6c57405SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "              RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);}
10897adad957SLisandro Dalcin   ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm,"              [%d]  %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr);
10907adad957SLisandro Dalcin   ierr = PetscSynchronizedFlush(((PetscObject)A)->comm);
1091f6c57405SHong Zhang 
1092f6c57405SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "              INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr);}
10937adad957SLisandro Dalcin   ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm,"              [%d] %d \n",lu->myid,lu->id.INFO(15));CHKERRQ(ierr);
10947adad957SLisandro Dalcin   ierr = PetscSynchronizedFlush(((PetscObject)A)->comm);
1095f6c57405SHong Zhang 
1096f6c57405SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "              INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr);}
10977adad957SLisandro Dalcin   ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm,"              [%d] %d \n",lu->myid,lu->id.INFO(16));CHKERRQ(ierr);
10987adad957SLisandro Dalcin   ierr = PetscSynchronizedFlush(((PetscObject)A)->comm);
1099f6c57405SHong Zhang 
1100f6c57405SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "              INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);}
11017adad957SLisandro Dalcin   ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm,"              [%d] %d \n",lu->myid,lu->id.INFO(23));CHKERRQ(ierr);
11027adad957SLisandro Dalcin   ierr = PetscSynchronizedFlush(((PetscObject)A)->comm);
1103f6c57405SHong Zhang 
1104f6c57405SHong Zhang   if (!lu->myid){ /* information from the host */
1105f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr);
1106f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr);
1107f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr);
1108f6c57405SHong Zhang 
1109f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr);
1110f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr);
1111f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr);
1112f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr);
1113f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr);
1114f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr);
1115f6c57405SHong 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);
1116f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",lu->id.INFOG(10));CHKERRQ(ierr);
1117f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",lu->id.INFOG(11));CHKERRQ(ierr);
1118f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr);
1119f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr);
1120f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr);
1121f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr);
1122f6c57405SHong 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);
1123f6c57405SHong 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);
1124f6c57405SHong 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);
1125f6c57405SHong 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);
1126f6c57405SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr);
1127f6c57405SHong 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);
1128f6c57405SHong 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);
1129f6c57405SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",lu->id.INFOG(23));CHKERRQ(ierr);
1130f6c57405SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",lu->id.INFOG(24));CHKERRQ(ierr);
1131f6c57405SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",lu->id.INFOG(25));CHKERRQ(ierr);
1132f6c57405SHong Zhang   }
1133f6c57405SHong Zhang   PetscFunctionReturn(0);
1134f6c57405SHong Zhang }
1135f6c57405SHong Zhang 
1136f6c57405SHong Zhang #undef __FUNCT__
1137f6c57405SHong Zhang #define __FUNCT__ "MatView_MUMPS"
1138b24902e0SBarry Smith PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
1139b24902e0SBarry Smith {
1140f6c57405SHong Zhang   PetscErrorCode    ierr;
1141f6c57405SHong Zhang   PetscTruth        iascii;
1142f6c57405SHong Zhang   PetscViewerFormat format;
1143cb828f0fSHong Zhang   Mat_MUMPS         *mumps=(Mat_MUMPS*)A->spptr;
1144f6c57405SHong Zhang 
1145f6c57405SHong Zhang   PetscFunctionBegin;
1146cb828f0fSHong Zhang   /* check if matrix is mumps type */
1147cb828f0fSHong Zhang   if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0);
1148cb828f0fSHong Zhang 
1149f6c57405SHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1150f6c57405SHong Zhang   if (iascii) {
1151f6c57405SHong Zhang     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1152cb828f0fSHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO || mumps->mumpsview){
1153cb828f0fSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
1154cb828f0fSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                  %d \n",mumps->id.sym);CHKERRQ(ierr);
1155cb828f0fSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):           %d \n",mumps->id.par);CHKERRQ(ierr);
1156cb828f0fSHong Zhang       if (mumps->mumpsview){ /* View all MUMPS parameters */
1157f6c57405SHong Zhang         ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr);
1158f6c57405SHong Zhang       }
1159f6c57405SHong Zhang     }
1160cb828f0fSHong Zhang   }
1161f6c57405SHong Zhang   PetscFunctionReturn(0);
1162f6c57405SHong Zhang }
1163f6c57405SHong Zhang 
116435bd34faSBarry Smith #undef __FUNCT__
116535bd34faSBarry Smith #define __FUNCT__ "MatGetInfo_MUMPS"
116635bd34faSBarry Smith PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
116735bd34faSBarry Smith {
1168cb828f0fSHong Zhang   Mat_MUMPS  *mumps =(Mat_MUMPS*)A->spptr;
116935bd34faSBarry Smith 
117035bd34faSBarry Smith   PetscFunctionBegin;
117135bd34faSBarry Smith   info->block_size        = 1.0;
1172cb828f0fSHong Zhang   info->nz_allocated      = mumps->id.INFOG(20);
1173cb828f0fSHong Zhang   info->nz_used           = mumps->id.INFOG(20);
117435bd34faSBarry Smith   info->nz_unneeded       = 0.0;
117535bd34faSBarry Smith   info->assemblies        = 0.0;
117635bd34faSBarry Smith   info->mallocs           = 0.0;
117735bd34faSBarry Smith   info->memory            = 0.0;
117835bd34faSBarry Smith   info->fill_ratio_given  = 0;
117935bd34faSBarry Smith   info->fill_ratio_needed = 0;
118035bd34faSBarry Smith   info->factor_mallocs    = 0;
118135bd34faSBarry Smith   PetscFunctionReturn(0);
118235bd34faSBarry Smith }
118335bd34faSBarry Smith 
118424b6179bSKris Buschelman /*MC
118541c8de11SBarry Smith   MAT_SOLVER_MUMPS -  A matrix type providing direct solvers (LU and Cholesky) for
118624b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
118724b6179bSKris Buschelman 
118841c8de11SBarry Smith   Works with MATAIJ and MATSBAIJ matrices
118924b6179bSKris Buschelman 
119024b6179bSKris Buschelman   Options Database Keys:
119141c8de11SBarry Smith + -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
119224b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,...,4> - print level
119324b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
119424b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
119524b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
119624b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
119794b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
119824b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
119924b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
120024b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
120124b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
120224b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
120324b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
120424b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
120524b6179bSKris Buschelman 
120624b6179bSKris Buschelman   Level: beginner
120724b6179bSKris Buschelman 
120841c8de11SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage
120941c8de11SBarry Smith 
121024b6179bSKris Buschelman M*/
121124b6179bSKris Buschelman 
12122877fffaSHong Zhang EXTERN_C_BEGIN
121335bd34faSBarry Smith #undef __FUNCT__
121435bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_mumps"
121535bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type)
121635bd34faSBarry Smith {
121735bd34faSBarry Smith   PetscFunctionBegin;
121835bd34faSBarry Smith   *type = MAT_SOLVER_MUMPS;
121935bd34faSBarry Smith   PetscFunctionReturn(0);
122035bd34faSBarry Smith }
122135bd34faSBarry Smith EXTERN_C_END
122235bd34faSBarry Smith 
122335bd34faSBarry Smith EXTERN_C_BEGIN
1224bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI AIJ matrices */
12252877fffaSHong Zhang #undef __FUNCT__
1226bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_aij_mumps"
1227bccb9932SShri Abhyankar PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
12282877fffaSHong Zhang {
12292877fffaSHong Zhang   Mat            B;
12302877fffaSHong Zhang   PetscErrorCode ierr;
12312877fffaSHong Zhang   Mat_MUMPS      *mumps;
1232bccb9932SShri Abhyankar   PetscTruth     isSeqAIJ;
12332877fffaSHong Zhang 
12342877fffaSHong Zhang   PetscFunctionBegin;
12352877fffaSHong Zhang   /* Create the factorization matrix */
1236bccb9932SShri Abhyankar   ierr = PetscTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr);
12372877fffaSHong Zhang   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
12382877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
12392877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1240bccb9932SShri Abhyankar   if (isSeqAIJ) {
12412877fffaSHong Zhang     ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
1242bccb9932SShri Abhyankar   } else {
1243bccb9932SShri Abhyankar     ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
1244bccb9932SShri Abhyankar   }
12452877fffaSHong Zhang 
124616ebf90aSShri Abhyankar   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
12472877fffaSHong Zhang   B->ops->view             = MatView_MUMPS;
124835bd34faSBarry Smith   B->ops->getinfo          = MatGetInfo_MUMPS;
124935bd34faSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
125097969023SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl",MatMumpsSetIcntl);CHKERRQ(ierr);
1251450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
1252450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
1253d5f3da31SBarry Smith     B->factortype = MAT_FACTOR_LU;
1254bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
1255bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
1256dcd589f8SShri Abhyankar   } else {
125767877ebaSShri Abhyankar     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
1258450b117fSShri Abhyankar     B->factortype = MAT_FACTOR_CHOLESKY;
1259bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
1260bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
1261450b117fSShri Abhyankar   }
12622877fffaSHong Zhang 
12632877fffaSHong Zhang   mumps->CleanUpMUMPS              = PETSC_FALSE;
12642877fffaSHong Zhang   mumps->isAIJ                     = PETSC_TRUE;
12652877fffaSHong Zhang   mumps->scat_rhs                  = PETSC_NULL;
12662877fffaSHong Zhang   mumps->scat_sol                  = PETSC_NULL;
12672877fffaSHong Zhang   mumps->nSolve                    = 0;
12682877fffaSHong Zhang   mumps->MatDestroy                = B->ops->destroy;
12692877fffaSHong Zhang   B->ops->destroy                  = MatDestroy_MUMPS;
12702877fffaSHong Zhang   B->spptr                         = (void*)mumps;
12712877fffaSHong Zhang 
12722877fffaSHong Zhang   *F = B;
12732877fffaSHong Zhang   PetscFunctionReturn(0);
12742877fffaSHong Zhang }
12752877fffaSHong Zhang EXTERN_C_END
12762877fffaSHong Zhang 
1277bccb9932SShri Abhyankar 
12782877fffaSHong Zhang EXTERN_C_BEGIN
1279bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI SBAIJ matrices */
12802877fffaSHong Zhang #undef __FUNCT__
1281bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_sbaij_mumps"
1282bccb9932SShri Abhyankar PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
12832877fffaSHong Zhang {
12842877fffaSHong Zhang   Mat            B;
12852877fffaSHong Zhang   PetscErrorCode ierr;
12862877fffaSHong Zhang   Mat_MUMPS      *mumps;
1287bccb9932SShri Abhyankar   PetscTruth     isSeqSBAIJ;
12882877fffaSHong Zhang 
12892877fffaSHong Zhang   PetscFunctionBegin;
1290e7e72b3dSBarry Smith   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
1291bccb9932SShri 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");
1292bccb9932SShri Abhyankar   ierr = PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
12932877fffaSHong Zhang   /* Create the factorization matrix */
12942877fffaSHong Zhang   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
12952877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
12962877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
129716ebf90aSShri Abhyankar   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
1298bccb9932SShri Abhyankar   if (isSeqSBAIJ) {
1299bccb9932SShri Abhyankar     ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
130016ebf90aSShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
1301dcd589f8SShri Abhyankar   } else {
1302bccb9932SShri Abhyankar     ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
1303bccb9932SShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
1304bccb9932SShri Abhyankar   }
1305bccb9932SShri Abhyankar 
130667877ebaSShri Abhyankar   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
1307bccb9932SShri Abhyankar   B->ops->view                   = MatView_MUMPS;
1308bccb9932SShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
1309bccb9932SShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl",MatMumpsSetIcntl);CHKERRQ(ierr);
1310f4762488SHong Zhang   B->factortype                   = MAT_FACTOR_CHOLESKY;
1311a214ac2aSShri Abhyankar 
1312a214ac2aSShri Abhyankar   mumps->CleanUpMUMPS              = PETSC_FALSE;
1313bccb9932SShri Abhyankar   mumps->isAIJ                     = PETSC_FALSE;
13142877fffaSHong Zhang   mumps->scat_rhs                  = PETSC_NULL;
13152877fffaSHong Zhang   mumps->scat_sol                  = PETSC_NULL;
13162877fffaSHong Zhang   mumps->nSolve                    = 0;
1317f3c0ef26SHong Zhang   mumps->MatDestroy                = B->ops->destroy;
1318f3c0ef26SHong Zhang   B->ops->destroy                  = MatDestroy_MUMPS;
13192877fffaSHong Zhang   B->spptr                         = (void*)mumps;
1320f3c0ef26SHong Zhang 
13212877fffaSHong Zhang   *F = B;
13222877fffaSHong Zhang   PetscFunctionReturn(0);
13232877fffaSHong Zhang }
13242877fffaSHong Zhang EXTERN_C_END
132597969023SHong Zhang 
1326450b117fSShri Abhyankar EXTERN_C_BEGIN
1327450b117fSShri Abhyankar #undef __FUNCT__
1328bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_baij_mumps"
1329bccb9932SShri Abhyankar PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
133067877ebaSShri Abhyankar {
133167877ebaSShri Abhyankar   Mat            B;
133267877ebaSShri Abhyankar   PetscErrorCode ierr;
133367877ebaSShri Abhyankar   Mat_MUMPS      *mumps;
1334bccb9932SShri Abhyankar   PetscTruth     isSeqBAIJ;
133567877ebaSShri Abhyankar 
133667877ebaSShri Abhyankar   PetscFunctionBegin;
133767877ebaSShri Abhyankar   /* Create the factorization matrix */
1338bccb9932SShri Abhyankar   ierr = PetscTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr);
133967877ebaSShri Abhyankar   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
134067877ebaSShri Abhyankar   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
134167877ebaSShri Abhyankar   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1342bccb9932SShri Abhyankar   if (isSeqBAIJ) {
134367877ebaSShri Abhyankar     ierr = MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,PETSC_NULL);CHKERRQ(ierr);
1344bccb9932SShri Abhyankar   } else {
134567877ebaSShri Abhyankar     ierr = MatMPIBAIJSetPreallocation(B,A->rmap->bs,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
1346bccb9932SShri Abhyankar   }
1347450b117fSShri Abhyankar 
134867877ebaSShri Abhyankar   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
1349450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
1350450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
1351450b117fSShri Abhyankar     B->factortype = MAT_FACTOR_LU;
1352bccb9932SShri Abhyankar     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
1353bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
1354450b117fSShri Abhyankar   }
1355bccb9932SShri Abhyankar   else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n");
1356bccb9932SShri Abhyankar 
1357450b117fSShri Abhyankar   B->ops->view             = MatView_MUMPS;
1358450b117fSShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
1359450b117fSShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl",MatMumpsSetIcntl);CHKERRQ(ierr);
1360450b117fSShri Abhyankar 
1361450b117fSShri Abhyankar   mumps->CleanUpMUMPS              = PETSC_FALSE;
1362450b117fSShri Abhyankar   mumps->isAIJ                     = PETSC_TRUE;
1363450b117fSShri Abhyankar   mumps->scat_rhs                  = PETSC_NULL;
1364450b117fSShri Abhyankar   mumps->scat_sol                  = PETSC_NULL;
1365450b117fSShri Abhyankar   mumps->nSolve                    = 0;
1366450b117fSShri Abhyankar   mumps->MatDestroy                = B->ops->destroy;
1367450b117fSShri Abhyankar   B->ops->destroy                  = MatDestroy_MUMPS;
1368450b117fSShri Abhyankar   B->spptr                         = (void*)mumps;
1369450b117fSShri Abhyankar 
1370450b117fSShri Abhyankar   *F = B;
1371450b117fSShri Abhyankar   PetscFunctionReturn(0);
1372450b117fSShri Abhyankar }
1373450b117fSShri Abhyankar EXTERN_C_END
1374a214ac2aSShri Abhyankar 
137561288eaaSHong Zhang /* -------------------------------------------------------------------------------------------*/
137661288eaaSHong Zhang /*@
137761288eaaSHong Zhang   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()
137861288eaaSHong Zhang 
137961288eaaSHong Zhang    Collective on Mat
138061288eaaSHong Zhang 
138161288eaaSHong Zhang    Input Parameters:
138261288eaaSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
138361288eaaSHong Zhang .  idx - index of MUMPS parameter array ICNTL()
138461288eaaSHong Zhang -  icntl - value of MUMPS ICNTL(imumps)
138561288eaaSHong Zhang 
138661288eaaSHong Zhang   Options Database:
138761288eaaSHong Zhang .   -mat_mumps_icntl_<idx> <icntl>
138861288eaaSHong Zhang 
138961288eaaSHong Zhang    Level: beginner
139061288eaaSHong Zhang 
139161288eaaSHong Zhang    References: MUMPS Users' Guide
139261288eaaSHong Zhang 
139361288eaaSHong Zhang .seealso: MatGetFactor()
139461288eaaSHong Zhang @*/
139597969023SHong Zhang #undef __FUNCT__
139697969023SHong Zhang #define __FUNCT__ "MatMumpsSetIcntl"
139786e5884dSMatthew Knepley PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt idx,PetscInt icntl)
139897969023SHong Zhang {
1399dcd589f8SShri Abhyankar   Mat_MUMPS      *lu =(Mat_MUMPS*)(F)->spptr;
140097969023SHong Zhang 
140197969023SHong Zhang   PetscFunctionBegin;
140261288eaaSHong Zhang   lu->id.ICNTL(idx) = icntl;
140397969023SHong Zhang   PetscFunctionReturn(0);
140497969023SHong Zhang }
140597969023SHong Zhang 
1406