xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision fb8376fbe8eb249cd5db38acdecbce5701dbf5d8)
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
213d472b54SHong Zhang #define JOB_FACTSYMBOLIC 1
223d472b54SHong Zhang #define JOB_FACTNUMERIC 2
233d472b54SHong Zhang #define JOB_SOLVE 3
24397b6df1SKris Buschelman #define JOB_END -2
253d472b54SHong Zhang 
263d472b54SHong Zhang 
27397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */
28397b6df1SKris Buschelman #define ICNTL(I) icntl[(I)-1]
29397b6df1SKris Buschelman #define CNTL(I) cntl[(I)-1]
30397b6df1SKris Buschelman #define INFOG(I) infog[(I)-1]
31a7aca84bSHong Zhang #define INFO(I) info[(I)-1]
32397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1]
33adc1d99fSHong Zhang #define RINFO(I) rinfo[(I)-1]
34397b6df1SKris Buschelman 
35397b6df1SKris Buschelman typedef struct {
36397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
37397b6df1SKris Buschelman   ZMUMPS_STRUC_C id;
38397b6df1SKris Buschelman #else
39397b6df1SKris Buschelman   DMUMPS_STRUC_C id;
40397b6df1SKris Buschelman #endif
41397b6df1SKris Buschelman   MatStructure   matstruc;
42c1490034SHong Zhang   PetscMPIInt    myid,size;
4316ebf90aSShri Abhyankar   PetscInt       *irn,*jcn,nz,sym,nSolve;
44397b6df1SKris Buschelman   PetscScalar    *val;
45397b6df1SKris Buschelman   MPI_Comm       comm_mumps;
46329ec9b3SHong Zhang   VecScatter     scat_rhs, scat_sol;
47cb828f0fSHong Zhang   PetscTruth     isAIJ,CleanUpMUMPS,mumpsview;
48329ec9b3SHong Zhang   Vec            b_seq,x_seq;
4967334b25SHong Zhang   PetscErrorCode (*MatDestroy)(Mat);
50bccb9932SShri Abhyankar   PetscErrorCode (*ConvertToTriples)(Mat, int, MatReuse, int*, int**, int**, PetscScalar**);
51f0c56d0fSKris Buschelman } Mat_MUMPS;
52f0c56d0fSKris Buschelman 
53dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);
54b24902e0SBarry Smith 
5567877ebaSShri Abhyankar 
5667877ebaSShri Abhyankar /* MatConvertToTriples_A_B */
5767877ebaSShri Abhyankar /*convert Petsc matrix to triples: row[nz], col[nz], val[nz] */
58397b6df1SKris Buschelman /*
59397b6df1SKris Buschelman   input:
6067877ebaSShri Abhyankar     A       - matrix in aij,baij or sbaij (bs=1) format
61397b6df1SKris Buschelman     shift   - 0: C style output triple; 1: Fortran style output triple.
62bccb9932SShri Abhyankar     reuse   - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple
63bccb9932SShri Abhyankar               MAT_REUSE_MATRIX:   only the values in v array are updated
64397b6df1SKris Buschelman   output:
65397b6df1SKris Buschelman     nnz     - dim of r, c, and v (number of local nonzero entries of A)
66397b6df1SKris Buschelman     r, c, v - row and col index, matrix values (matrix triples)
67397b6df1SKris Buschelman  */
6816ebf90aSShri Abhyankar 
6916ebf90aSShri Abhyankar #undef __FUNCT__
7016ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqaij"
71bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
72b24902e0SBarry Smith {
73185f6596SHong Zhang   const PetscInt   *ai,*aj,*ajj,M=A->rmap->n;
7467877ebaSShri Abhyankar   PetscInt         nz,rnz,i,j;
75dfbe8321SBarry Smith   PetscErrorCode   ierr;
76c1490034SHong Zhang   PetscInt         *row,*col;
7716ebf90aSShri Abhyankar   Mat_SeqAIJ       *aa=(Mat_SeqAIJ*)A->data;
78397b6df1SKris Buschelman 
79397b6df1SKris Buschelman   PetscFunctionBegin;
8016ebf90aSShri Abhyankar   *v=aa->a;
81bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX){
8216ebf90aSShri Abhyankar     nz = aa->nz; ai = aa->i; aj = aa->j;
8316ebf90aSShri Abhyankar     *nnz = nz;
84185f6596SHong Zhang     ierr = PetscMalloc(2*nz*sizeof(PetscInt), &row);CHKERRQ(ierr);
85185f6596SHong Zhang     col  = row + nz;
86185f6596SHong Zhang 
8716ebf90aSShri Abhyankar     nz = 0;
8816ebf90aSShri Abhyankar     for(i=0; i<M; i++) {
8916ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
9067877ebaSShri Abhyankar       ajj = aj + ai[i];
9167877ebaSShri Abhyankar       for(j=0; j<rnz; j++) {
9267877ebaSShri Abhyankar 	row[nz] = i+shift; col[nz++] = ajj[j] + shift;
9316ebf90aSShri Abhyankar       }
9416ebf90aSShri Abhyankar     }
9516ebf90aSShri Abhyankar     *r = row; *c = col;
9616ebf90aSShri Abhyankar   }
9716ebf90aSShri Abhyankar   PetscFunctionReturn(0);
9816ebf90aSShri Abhyankar }
99397b6df1SKris Buschelman 
10016ebf90aSShri Abhyankar #undef __FUNCT__
10167877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqbaij_seqaij"
102bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
10367877ebaSShri Abhyankar {
10467877ebaSShri Abhyankar   Mat_SeqBAIJ        *aa=(Mat_SeqBAIJ*)A->data;
10567877ebaSShri Abhyankar   const PetscInt     *ai,*aj,*ajj,bs=A->rmap->bs,bs2=aa->bs2,M=A->rmap->N/bs;
10667877ebaSShri Abhyankar   PetscInt           nz,idx=0,rnz,i,j,k,m,ii;
10767877ebaSShri Abhyankar   PetscErrorCode     ierr;
10867877ebaSShri Abhyankar   PetscInt           *row,*col;
10967877ebaSShri Abhyankar 
11067877ebaSShri Abhyankar   PetscFunctionBegin;
111cf3759fdSShri Abhyankar   *v = aa->a;
112bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX){
113cf3759fdSShri Abhyankar     ai = aa->i; aj = aa->j;
11467877ebaSShri Abhyankar     nz = bs2*aa->nz;
11567877ebaSShri Abhyankar     *nnz = nz;
116185f6596SHong Zhang     ierr = PetscMalloc(2*nz*sizeof(PetscInt), &row);CHKERRQ(ierr);
117185f6596SHong Zhang     col  = row + nz;
118185f6596SHong Zhang 
11967877ebaSShri Abhyankar     for(i=0; i<M; i++) {
12067877ebaSShri Abhyankar       ii = 0;
12167877ebaSShri Abhyankar       ajj = aj + ai[i];
12267877ebaSShri Abhyankar       rnz = ai[i+1] - ai[i];
12367877ebaSShri Abhyankar       for(k=0; k<rnz; k++) {
12467877ebaSShri Abhyankar 	for(j=0; j<bs; j++) {
12567877ebaSShri Abhyankar 	  for(m=0; m<bs; m++) {
12667877ebaSShri Abhyankar 	    row[idx]     = i*bs + m + shift;
127cf3759fdSShri Abhyankar 	    col[idx++]   = bs*(ajj[k]) + j + shift;
12867877ebaSShri Abhyankar 	  }
12967877ebaSShri Abhyankar 	}
13067877ebaSShri Abhyankar       }
13167877ebaSShri Abhyankar     }
132cf3759fdSShri Abhyankar     *r = row; *c = col;
13367877ebaSShri Abhyankar   }
13467877ebaSShri Abhyankar   PetscFunctionReturn(0);
13567877ebaSShri Abhyankar }
13667877ebaSShri Abhyankar 
13767877ebaSShri Abhyankar #undef __FUNCT__
13816ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqsbaij_seqsbaij"
139bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
14016ebf90aSShri Abhyankar {
14167877ebaSShri Abhyankar   const PetscInt   *ai, *aj,*ajj,M=A->rmap->n;
14267877ebaSShri Abhyankar   PetscInt         nz,rnz,i,j;
14316ebf90aSShri Abhyankar   PetscErrorCode   ierr;
14416ebf90aSShri Abhyankar   PetscInt         *row,*col;
14516ebf90aSShri Abhyankar   Mat_SeqSBAIJ     *aa=(Mat_SeqSBAIJ*)A->data;
14616ebf90aSShri Abhyankar 
14716ebf90aSShri Abhyankar   PetscFunctionBegin;
148bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX){
14916ebf90aSShri Abhyankar     nz = aa->nz;ai=aa->i; aj=aa->j;*v=aa->a;
15016ebf90aSShri Abhyankar     *nnz = nz;
151185f6596SHong Zhang     ierr = PetscMalloc(2*nz*sizeof(PetscInt), &row);CHKERRQ(ierr);
152185f6596SHong Zhang     col  = row + nz;
153185f6596SHong Zhang 
15416ebf90aSShri Abhyankar     nz = 0;
15516ebf90aSShri Abhyankar     for(i=0; i<M; i++) {
15616ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
15767877ebaSShri Abhyankar       ajj = aj + ai[i];
15867877ebaSShri Abhyankar       for(j=0; j<rnz; j++) {
15967877ebaSShri Abhyankar 	row[nz] = i+shift; col[nz++] = ajj[j] + shift;
16016ebf90aSShri Abhyankar       }
16116ebf90aSShri Abhyankar     }
16216ebf90aSShri Abhyankar     *r = row; *c = col;
16316ebf90aSShri Abhyankar   }
16416ebf90aSShri Abhyankar   PetscFunctionReturn(0);
16516ebf90aSShri Abhyankar }
16616ebf90aSShri Abhyankar 
16716ebf90aSShri Abhyankar #undef __FUNCT__
16816ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqsbaij"
169bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
17016ebf90aSShri Abhyankar {
17167877ebaSShri Abhyankar   const PetscInt     *ai,*aj,*ajj,*adiag,M=A->rmap->n;
17267877ebaSShri Abhyankar   PetscInt           nz,rnz,i,j;
17367877ebaSShri Abhyankar   const PetscScalar  *av,*v1;
17416ebf90aSShri Abhyankar   PetscScalar        *val;
17516ebf90aSShri Abhyankar   PetscErrorCode     ierr;
17616ebf90aSShri Abhyankar   PetscInt           *row,*col;
17716ebf90aSShri Abhyankar   Mat_SeqSBAIJ       *aa=(Mat_SeqSBAIJ*)A->data;
17816ebf90aSShri Abhyankar 
17916ebf90aSShri Abhyankar   PetscFunctionBegin;
18016ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j;av=aa->a;
18116ebf90aSShri Abhyankar   adiag=aa->diag;
182bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX){
18316ebf90aSShri Abhyankar     nz = M + (aa->nz-M)/2;
18416ebf90aSShri Abhyankar     *nnz = nz;
185185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
186185f6596SHong Zhang     col  = row + nz;
187185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
188185f6596SHong Zhang 
18916ebf90aSShri Abhyankar     nz = 0;
19016ebf90aSShri Abhyankar     for(i=0; i<M; i++) {
19116ebf90aSShri Abhyankar       rnz = ai[i+1] - adiag[i];
19267877ebaSShri Abhyankar       ajj  = aj + adiag[i];
193cf3759fdSShri Abhyankar       v1   = av + adiag[i];
19467877ebaSShri Abhyankar       for(j=0; j<rnz; j++) {
19567877ebaSShri Abhyankar 	row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
19616ebf90aSShri Abhyankar       }
19716ebf90aSShri Abhyankar     }
19816ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
199397b6df1SKris Buschelman   } else {
20016ebf90aSShri Abhyankar     nz = 0; val = *v;
20116ebf90aSShri Abhyankar     for(i=0; i <M; i++) {
20216ebf90aSShri Abhyankar       rnz = ai[i+1] - adiag[i];
20367877ebaSShri Abhyankar       ajj = aj + adiag[i];
20467877ebaSShri Abhyankar       v1  = av + adiag[i];
20567877ebaSShri Abhyankar       for(j=0; j<rnz; j++) {
20667877ebaSShri Abhyankar 	val[nz++] = v1[j];
20716ebf90aSShri Abhyankar       }
20816ebf90aSShri Abhyankar     }
20916ebf90aSShri Abhyankar   }
21016ebf90aSShri Abhyankar   PetscFunctionReturn(0);
21116ebf90aSShri Abhyankar }
21216ebf90aSShri Abhyankar 
21316ebf90aSShri Abhyankar #undef __FUNCT__
21416ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpisbaij_mpisbaij"
215bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
21616ebf90aSShri Abhyankar {
21716ebf90aSShri Abhyankar   const PetscInt     *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
21816ebf90aSShri Abhyankar   PetscErrorCode     ierr;
21916ebf90aSShri Abhyankar   PetscInt           rstart,nz,i,j,jj,irow,countA,countB;
22016ebf90aSShri Abhyankar   PetscInt           *row,*col;
22116ebf90aSShri Abhyankar   const PetscScalar  *av, *bv,*v1,*v2;
22216ebf90aSShri Abhyankar   PetscScalar        *val;
223397b6df1SKris Buschelman   Mat_MPISBAIJ       *mat =  (Mat_MPISBAIJ*)A->data;
224397b6df1SKris Buschelman   Mat_SeqSBAIJ       *aa=(Mat_SeqSBAIJ*)(mat->A)->data;
225397b6df1SKris Buschelman   Mat_SeqBAIJ        *bb=(Mat_SeqBAIJ*)(mat->B)->data;
22616ebf90aSShri Abhyankar 
22716ebf90aSShri Abhyankar   PetscFunctionBegin;
228d0f46423SBarry Smith   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
229397b6df1SKris Buschelman   garray = mat->garray;
230397b6df1SKris Buschelman   av=aa->a; bv=bb->a;
231397b6df1SKris Buschelman 
232bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX){
23316ebf90aSShri Abhyankar     nz = aa->nz + bb->nz;
23416ebf90aSShri Abhyankar     *nnz = nz;
235185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
236185f6596SHong Zhang     col  = row + nz;
237185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
238185f6596SHong Zhang 
239397b6df1SKris Buschelman     *r = row; *c = col; *v = val;
240397b6df1SKris Buschelman   } else {
241397b6df1SKris Buschelman     row = *r; col = *c; val = *v;
242397b6df1SKris Buschelman   }
243397b6df1SKris Buschelman 
244028e57e8SHong Zhang   jj = 0; irow = rstart;
245397b6df1SKris Buschelman   for ( i=0; i<m; i++ ) {
246397b6df1SKris Buschelman     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
247397b6df1SKris Buschelman     countA = ai[i+1] - ai[i];
248397b6df1SKris Buschelman     countB = bi[i+1] - bi[i];
249397b6df1SKris Buschelman     bjj    = bj + bi[i];
25016ebf90aSShri Abhyankar     v1     = av + ai[i];
25116ebf90aSShri Abhyankar     v2     = bv + bi[i];
252397b6df1SKris Buschelman 
253397b6df1SKris Buschelman     /* A-part */
254397b6df1SKris Buschelman     for (j=0; j<countA; j++){
255bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
256397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
257397b6df1SKris Buschelman       }
25816ebf90aSShri Abhyankar       val[jj++] = v1[j];
259397b6df1SKris Buschelman     }
26016ebf90aSShri Abhyankar 
26116ebf90aSShri Abhyankar     /* B-part */
26216ebf90aSShri Abhyankar     for(j=0; j < countB; j++){
263bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
264397b6df1SKris Buschelman 	row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
265397b6df1SKris Buschelman       }
26616ebf90aSShri Abhyankar       val[jj++] = v2[j];
26716ebf90aSShri Abhyankar     }
26816ebf90aSShri Abhyankar     irow++;
26916ebf90aSShri Abhyankar   }
27016ebf90aSShri Abhyankar   PetscFunctionReturn(0);
27116ebf90aSShri Abhyankar }
27216ebf90aSShri Abhyankar 
27316ebf90aSShri Abhyankar #undef __FUNCT__
27416ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpiaij"
275bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
27616ebf90aSShri Abhyankar {
27716ebf90aSShri Abhyankar   const PetscInt     *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
27816ebf90aSShri Abhyankar   PetscErrorCode     ierr;
27916ebf90aSShri Abhyankar   PetscInt           rstart,nz,i,j,jj,irow,countA,countB;
28016ebf90aSShri Abhyankar   PetscInt           *row,*col;
28116ebf90aSShri Abhyankar   const PetscScalar  *av, *bv,*v1,*v2;
28216ebf90aSShri Abhyankar   PetscScalar        *val;
28316ebf90aSShri Abhyankar   Mat_MPIAIJ         *mat =  (Mat_MPIAIJ*)A->data;
28416ebf90aSShri Abhyankar   Mat_SeqAIJ         *aa=(Mat_SeqAIJ*)(mat->A)->data;
28516ebf90aSShri Abhyankar   Mat_SeqAIJ         *bb=(Mat_SeqAIJ*)(mat->B)->data;
28616ebf90aSShri Abhyankar 
28716ebf90aSShri Abhyankar   PetscFunctionBegin;
28816ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
28916ebf90aSShri Abhyankar   garray = mat->garray;
29016ebf90aSShri Abhyankar   av=aa->a; bv=bb->a;
29116ebf90aSShri Abhyankar 
292bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX){
29316ebf90aSShri Abhyankar     nz = aa->nz + bb->nz;
29416ebf90aSShri Abhyankar     *nnz = nz;
295185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
296185f6596SHong Zhang     col  = row + nz;
297185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
298185f6596SHong Zhang 
29916ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
30016ebf90aSShri Abhyankar   } else {
30116ebf90aSShri Abhyankar     row = *r; col = *c; val = *v;
30216ebf90aSShri Abhyankar   }
30316ebf90aSShri Abhyankar 
30416ebf90aSShri Abhyankar   jj = 0; irow = rstart;
30516ebf90aSShri Abhyankar   for ( i=0; i<m; i++ ) {
30616ebf90aSShri Abhyankar     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
30716ebf90aSShri Abhyankar     countA = ai[i+1] - ai[i];
30816ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
30916ebf90aSShri Abhyankar     bjj    = bj + bi[i];
31016ebf90aSShri Abhyankar     v1     = av + ai[i];
31116ebf90aSShri Abhyankar     v2     = bv + bi[i];
31216ebf90aSShri Abhyankar 
31316ebf90aSShri Abhyankar     /* A-part */
31416ebf90aSShri Abhyankar     for (j=0; j<countA; j++){
315bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX){
31616ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
31716ebf90aSShri Abhyankar       }
31816ebf90aSShri Abhyankar       val[jj++] = v1[j];
31916ebf90aSShri Abhyankar     }
32016ebf90aSShri Abhyankar 
32116ebf90aSShri Abhyankar     /* B-part */
32216ebf90aSShri Abhyankar     for(j=0; j < countB; j++){
323bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX){
32416ebf90aSShri Abhyankar 	row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
32516ebf90aSShri Abhyankar       }
32616ebf90aSShri Abhyankar       val[jj++] = v2[j];
32716ebf90aSShri Abhyankar     }
32816ebf90aSShri Abhyankar     irow++;
32916ebf90aSShri Abhyankar   }
33016ebf90aSShri Abhyankar   PetscFunctionReturn(0);
33116ebf90aSShri Abhyankar }
33216ebf90aSShri Abhyankar 
33316ebf90aSShri Abhyankar #undef __FUNCT__
33467877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpibaij_mpiaij"
335bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
33667877ebaSShri Abhyankar {
33767877ebaSShri Abhyankar   Mat_MPIBAIJ        *mat =  (Mat_MPIBAIJ*)A->data;
33867877ebaSShri Abhyankar   Mat_SeqBAIJ        *aa=(Mat_SeqBAIJ*)(mat->A)->data;
33967877ebaSShri Abhyankar   Mat_SeqBAIJ        *bb=(Mat_SeqBAIJ*)(mat->B)->data;
34067877ebaSShri Abhyankar   const PetscInt     *ai = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj;
34167877ebaSShri Abhyankar   const PetscInt     *garray = mat->garray,mbs=mat->mbs,rstartbs=mat->rstartbs;
34267877ebaSShri Abhyankar   const PetscInt     bs = A->rmap->bs,bs2=mat->bs2;
34367877ebaSShri Abhyankar   PetscErrorCode     ierr;
34467877ebaSShri Abhyankar   PetscInt           nz,i,j,k,n,jj,irow,countA,countB,idx;
34567877ebaSShri Abhyankar   PetscInt           *row,*col;
34667877ebaSShri Abhyankar   const PetscScalar  *av=aa->a, *bv=bb->a,*v1,*v2;
34767877ebaSShri Abhyankar   PetscScalar        *val;
34867877ebaSShri Abhyankar 
34967877ebaSShri Abhyankar   PetscFunctionBegin;
35067877ebaSShri Abhyankar 
351bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
35267877ebaSShri Abhyankar     nz = bs2*(aa->nz + bb->nz);
35367877ebaSShri Abhyankar     *nnz = nz;
354185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
355185f6596SHong Zhang     col  = row + nz;
356185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
357185f6596SHong Zhang 
35867877ebaSShri Abhyankar     *r = row; *c = col; *v = val;
35967877ebaSShri Abhyankar   } else {
36067877ebaSShri Abhyankar     row = *r; col = *c; val = *v;
36167877ebaSShri Abhyankar   }
36267877ebaSShri Abhyankar 
36367877ebaSShri Abhyankar   jj = 0; irow = rstartbs;
36467877ebaSShri Abhyankar   for ( i=0; i<mbs; i++ ) {
36567877ebaSShri Abhyankar     countA = ai[i+1] - ai[i];
36667877ebaSShri Abhyankar     countB = bi[i+1] - bi[i];
36767877ebaSShri Abhyankar     ajj    = aj + ai[i];
36867877ebaSShri Abhyankar     bjj    = bj + bi[i];
36967877ebaSShri Abhyankar     v1     = av + bs2*ai[i];
37067877ebaSShri Abhyankar     v2     = bv + bs2*bi[i];
37167877ebaSShri Abhyankar 
37267877ebaSShri Abhyankar     idx = 0;
37367877ebaSShri Abhyankar     /* A-part */
37467877ebaSShri Abhyankar     for (k=0; k<countA; k++){
37567877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
37667877ebaSShri Abhyankar 	for (n=0; n<bs; n++) {
377bccb9932SShri Abhyankar 	  if (reuse == MAT_INITIAL_MATRIX){
37867877ebaSShri Abhyankar 	    row[jj] = bs*irow + n + shift;
37967877ebaSShri Abhyankar 	    col[jj] = bs*(rstartbs + ajj[k]) + j + shift;
38067877ebaSShri Abhyankar 	  }
38167877ebaSShri Abhyankar 	  val[jj++] = v1[idx++];
38267877ebaSShri Abhyankar 	}
38367877ebaSShri Abhyankar       }
38467877ebaSShri Abhyankar     }
38567877ebaSShri Abhyankar 
38667877ebaSShri Abhyankar     idx = 0;
38767877ebaSShri Abhyankar     /* B-part */
38867877ebaSShri Abhyankar     for(k=0; k<countB; k++){
38967877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
39067877ebaSShri Abhyankar 	for (n=0; n<bs; n++) {
391bccb9932SShri Abhyankar 	  if (reuse == MAT_INITIAL_MATRIX){
39267877ebaSShri Abhyankar 	    row[jj] = bs*irow + n + shift;
39367877ebaSShri Abhyankar 	    col[jj] = bs*(garray[bjj[k]]) + j + shift;
39467877ebaSShri Abhyankar 	  }
39567877ebaSShri Abhyankar 	  val[jj++] = bv[idx++];
39667877ebaSShri Abhyankar 	}
39767877ebaSShri Abhyankar       }
39867877ebaSShri Abhyankar     }
39967877ebaSShri Abhyankar     irow++;
40067877ebaSShri Abhyankar   }
40167877ebaSShri Abhyankar   PetscFunctionReturn(0);
40267877ebaSShri Abhyankar }
40367877ebaSShri Abhyankar 
40467877ebaSShri Abhyankar #undef __FUNCT__
40516ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpisbaij"
406bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
40716ebf90aSShri Abhyankar {
40816ebf90aSShri Abhyankar   const PetscInt     *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
40916ebf90aSShri Abhyankar   PetscErrorCode     ierr;
41016ebf90aSShri Abhyankar   PetscInt           rstart,nz,nza,nzb_low,i,j,jj,irow,countA,countB;
41116ebf90aSShri Abhyankar   PetscInt           *row,*col;
41216ebf90aSShri Abhyankar   const PetscScalar  *av, *bv,*v1,*v2;
41316ebf90aSShri Abhyankar   PetscScalar        *val;
41416ebf90aSShri Abhyankar   Mat_MPIAIJ         *mat =  (Mat_MPIAIJ*)A->data;
41516ebf90aSShri Abhyankar   Mat_SeqAIJ         *aa=(Mat_SeqAIJ*)(mat->A)->data;
41616ebf90aSShri Abhyankar   Mat_SeqAIJ         *bb=(Mat_SeqAIJ*)(mat->B)->data;
41716ebf90aSShri Abhyankar 
41816ebf90aSShri Abhyankar   PetscFunctionBegin;
41916ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j; adiag=aa->diag;
42016ebf90aSShri Abhyankar   bi=bb->i; bj=bb->j; garray = mat->garray;
42116ebf90aSShri Abhyankar   av=aa->a; bv=bb->a;
42216ebf90aSShri Abhyankar   rstart = A->rmap->rstart;
42316ebf90aSShri Abhyankar 
424bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
42516ebf90aSShri Abhyankar     nza = 0;nzb_low = 0;
42616ebf90aSShri Abhyankar     for(i=0; i<m; i++){
42716ebf90aSShri Abhyankar       nza     = nza + (ai[i+1] - adiag[i]);
42816ebf90aSShri Abhyankar       countB  = bi[i+1] - bi[i];
42916ebf90aSShri Abhyankar       bjj     = bj + bi[i];
43016ebf90aSShri Abhyankar 
43116ebf90aSShri Abhyankar       j = 0;
43216ebf90aSShri Abhyankar       while(garray[bjj[j]] < rstart) {
43316ebf90aSShri Abhyankar 	if(j == countB) break;
43416ebf90aSShri Abhyankar 	j++;nzb_low++;
43516ebf90aSShri Abhyankar       }
43616ebf90aSShri Abhyankar     }
43716ebf90aSShri Abhyankar     /* Total nz = nz for the upper triangular A part + nz for the 2nd B part */
43816ebf90aSShri Abhyankar     nz = nza + (bb->nz - nzb_low);
43916ebf90aSShri Abhyankar     *nnz = nz;
440185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
441185f6596SHong Zhang     col  = row + nz;
442185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
443185f6596SHong Zhang 
44416ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
44516ebf90aSShri Abhyankar   } else {
44616ebf90aSShri Abhyankar     row = *r; col = *c; val = *v;
44716ebf90aSShri Abhyankar   }
44816ebf90aSShri Abhyankar 
44916ebf90aSShri Abhyankar   jj = 0; irow = rstart;
45016ebf90aSShri Abhyankar   for ( i=0; i<m; i++ ) {
45116ebf90aSShri Abhyankar     ajj    = aj + adiag[i];                 /* ptr to the beginning of the diagonal of this row */
45216ebf90aSShri Abhyankar     v1     = av + adiag[i];
45316ebf90aSShri Abhyankar     countA = ai[i+1] - adiag[i];
45416ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
45516ebf90aSShri Abhyankar     bjj    = bj + bi[i];
45616ebf90aSShri Abhyankar     v2     = bv + bi[i];
45716ebf90aSShri Abhyankar 
45816ebf90aSShri Abhyankar      /* A-part */
45916ebf90aSShri Abhyankar     for (j=0; j<countA; j++){
460bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
46116ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
46216ebf90aSShri Abhyankar       }
46316ebf90aSShri Abhyankar       val[jj++] = v1[j];
46416ebf90aSShri Abhyankar     }
46516ebf90aSShri Abhyankar 
46616ebf90aSShri Abhyankar     /* B-part */
46716ebf90aSShri Abhyankar     for(j=0; j < countB; j++){
46816ebf90aSShri Abhyankar       if (garray[bjj[j]] > rstart) {
469bccb9932SShri Abhyankar 	if (reuse == MAT_INITIAL_MATRIX) {
47016ebf90aSShri Abhyankar 	  row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
47116ebf90aSShri Abhyankar 	}
47216ebf90aSShri Abhyankar 	val[jj++] = v2[j];
47316ebf90aSShri Abhyankar       }
474397b6df1SKris Buschelman     }
475397b6df1SKris Buschelman     irow++;
476397b6df1SKris Buschelman   }
477397b6df1SKris Buschelman   PetscFunctionReturn(0);
478397b6df1SKris Buschelman }
479397b6df1SKris Buschelman 
480397b6df1SKris Buschelman #undef __FUNCT__
4813924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS"
482dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A)
483dfbe8321SBarry Smith {
484f0c56d0fSKris Buschelman   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
485dfbe8321SBarry Smith   PetscErrorCode ierr;
486c1490034SHong Zhang   PetscMPIInt    size=lu->size;
487b24902e0SBarry Smith 
488397b6df1SKris Buschelman   PetscFunctionBegin;
489397b6df1SKris Buschelman   if (lu->CleanUpMUMPS) {
490397b6df1SKris Buschelman     /* Terminate instance, deallocate memories */
491329ec9b3SHong Zhang     if (size > 1){
49268653410SBarry Smith       ierr = PetscFree2(lu->id.sol_loc,lu->id.isol_loc);CHKERRQ(ierr);
493329ec9b3SHong Zhang       ierr = VecScatterDestroy(lu->scat_rhs);CHKERRQ(ierr);
494329ec9b3SHong Zhang       ierr = VecDestroy(lu->b_seq);CHKERRQ(ierr);
4952750af12SHong Zhang       if (lu->nSolve && lu->scat_sol){ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);}
4962750af12SHong Zhang       if (lu->nSolve && lu->x_seq){ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);}
497450b117fSShri Abhyankar     }
498185f6596SHong Zhang     ierr = PetscFree(lu->irn);CHKERRQ(ierr);
499397b6df1SKris Buschelman     lu->id.job=JOB_END;
500397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
501397b6df1SKris Buschelman     zmumps_c(&lu->id);
502397b6df1SKris Buschelman #else
503397b6df1SKris Buschelman     dmumps_c(&lu->id);
504397b6df1SKris Buschelman #endif
505397b6df1SKris Buschelman     ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr);
506397b6df1SKris Buschelman   }
50797969023SHong Zhang   /* clear composed functions */
50897969023SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatFactorGetSolverPackage_C","",PETSC_NULL);CHKERRQ(ierr);
509a1f19f5aSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSetMumpsIcntl_C","",PETSC_NULL);CHKERRQ(ierr);
51067334b25SHong Zhang   ierr = (lu->MatDestroy)(A);CHKERRQ(ierr);
511397b6df1SKris Buschelman   PetscFunctionReturn(0);
512397b6df1SKris Buschelman }
513397b6df1SKris Buschelman 
514397b6df1SKris Buschelman #undef __FUNCT__
515f6c57405SHong Zhang #define __FUNCT__ "MatSolve_MUMPS"
516b24902e0SBarry Smith PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x)
517b24902e0SBarry Smith {
518f0c56d0fSKris Buschelman   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
519d54de34fSKris Buschelman   PetscScalar    *array;
52067877ebaSShri Abhyankar   Vec            b_seq;
521329ec9b3SHong Zhang   IS             is_iden,is_petsc;
522dfbe8321SBarry Smith   PetscErrorCode ierr;
523329ec9b3SHong Zhang   PetscInt       i;
524397b6df1SKris Buschelman 
525397b6df1SKris Buschelman   PetscFunctionBegin;
526329ec9b3SHong Zhang   lu->id.nrhs = 1;
52767877ebaSShri Abhyankar   b_seq = lu->b_seq;
528397b6df1SKris Buschelman   if (lu->size > 1){
529329ec9b3SHong Zhang     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
53067877ebaSShri Abhyankar     ierr = VecScatterBegin(lu->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
53167877ebaSShri Abhyankar     ierr = VecScatterEnd(lu->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
53267877ebaSShri Abhyankar     if (!lu->myid) {ierr = VecGetArray(b_seq,&array);CHKERRQ(ierr);}
533397b6df1SKris Buschelman   } else {  /* size == 1 */
534397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
535397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
536397b6df1SKris Buschelman   }
537397b6df1SKris Buschelman   if (!lu->myid) { /* define rhs on the host */
5388278f211SHong Zhang     lu->id.nrhs = 1;
539397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
540397b6df1SKris Buschelman     lu->id.rhs = (mumps_double_complex*)array;
541397b6df1SKris Buschelman #else
542397b6df1SKris Buschelman     lu->id.rhs = array;
543397b6df1SKris Buschelman #endif
544397b6df1SKris Buschelman   }
545397b6df1SKris Buschelman 
546397b6df1SKris Buschelman   /* solve phase */
547329ec9b3SHong Zhang   /*-------------*/
5483d472b54SHong Zhang   lu->id.job = JOB_SOLVE;
549397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
550397b6df1SKris Buschelman   zmumps_c(&lu->id);
551397b6df1SKris Buschelman #else
552397b6df1SKris Buschelman   dmumps_c(&lu->id);
553397b6df1SKris Buschelman #endif
55465e19b50SBarry 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));
555397b6df1SKris Buschelman 
556329ec9b3SHong Zhang   if (lu->size > 1) { /* convert mumps distributed solution to petsc mpi x */
557329ec9b3SHong Zhang     if (!lu->nSolve){ /* create scatter scat_sol */
558329ec9b3SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,lu->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */
559329ec9b3SHong Zhang       for (i=0; i<lu->id.lsol_loc; i++){
560329ec9b3SHong Zhang         lu->id.isol_loc[i] -= 1; /* change Fortran style to C style */
561397b6df1SKris Buschelman       }
562329ec9b3SHong Zhang       ierr = ISCreateGeneral(PETSC_COMM_SELF,lu->id.lsol_loc,lu->id.isol_loc,&is_petsc);CHKERRQ(ierr);  /* to */
563329ec9b3SHong Zhang       ierr = VecScatterCreate(lu->x_seq,is_iden,x,is_petsc,&lu->scat_sol);CHKERRQ(ierr);
564329ec9b3SHong Zhang       ierr = ISDestroy(is_iden);CHKERRQ(ierr);
565329ec9b3SHong Zhang       ierr = ISDestroy(is_petsc);CHKERRQ(ierr);
566397b6df1SKris Buschelman     }
567ca9f406cSSatish Balay     ierr = VecScatterBegin(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
568ca9f406cSSatish Balay     ierr = VecScatterEnd(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
569329ec9b3SHong Zhang   }
570329ec9b3SHong Zhang   lu->nSolve++;
571397b6df1SKris Buschelman   PetscFunctionReturn(0);
572397b6df1SKris Buschelman }
573397b6df1SKris Buschelman 
574ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX)
575a58c3f20SHong Zhang /*
576a58c3f20SHong Zhang   input:
577a58c3f20SHong Zhang    F:        numeric factor
578a58c3f20SHong Zhang   output:
579a58c3f20SHong Zhang    nneg:     total number of negative pivots
580a58c3f20SHong Zhang    nzero:    0
581a58c3f20SHong Zhang    npos:     (global dimension of F) - nneg
582a58c3f20SHong Zhang */
583a58c3f20SHong Zhang 
584a58c3f20SHong Zhang #undef __FUNCT__
585a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
586dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
587a58c3f20SHong Zhang {
588a58c3f20SHong Zhang   Mat_MUMPS      *lu =(Mat_MUMPS*)F->spptr;
589dfbe8321SBarry Smith   PetscErrorCode ierr;
590c1490034SHong Zhang   PetscMPIInt    size;
591a58c3f20SHong Zhang 
592a58c3f20SHong Zhang   PetscFunctionBegin;
5937adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)F)->comm,&size);CHKERRQ(ierr);
594bcb30aebSHong 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 */
59565e19b50SBarry 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));
596a58c3f20SHong Zhang   if (nneg){
597a58c3f20SHong Zhang     if (!lu->myid){
598a58c3f20SHong Zhang       *nneg = lu->id.INFOG(12);
599a58c3f20SHong Zhang     }
600bcb30aebSHong Zhang     ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr);
601a58c3f20SHong Zhang   }
602a58c3f20SHong Zhang   if (nzero) *nzero = 0;
603d0f46423SBarry Smith   if (npos)  *npos  = F->rmap->N - (*nneg);
604a58c3f20SHong Zhang   PetscFunctionReturn(0);
605a58c3f20SHong Zhang }
606ace3df97SHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */
607a58c3f20SHong Zhang 
608397b6df1SKris Buschelman #undef __FUNCT__
609f6c57405SHong Zhang #define __FUNCT__ "MatFactorNumeric_MUMPS"
6100481f469SBarry Smith PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
611af281ebdSHong Zhang {
612dcd589f8SShri Abhyankar   Mat_MUMPS       *lu =(Mat_MUMPS*)(F)->spptr;
6136849ba73SBarry Smith   PetscErrorCode  ierr;
614bccb9932SShri Abhyankar   MatReuse        reuse;
615e09efc27SHong Zhang   Mat             F_diag;
61616ebf90aSShri Abhyankar   PetscTruth      isMPIAIJ;
617397b6df1SKris Buschelman 
618397b6df1SKris Buschelman   PetscFunctionBegin;
619bccb9932SShri Abhyankar   reuse = MAT_REUSE_MATRIX;
620bccb9932SShri Abhyankar   ierr = (*lu->ConvertToTriples)(A, 1, reuse, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
621397b6df1SKris Buschelman 
622397b6df1SKris Buschelman   /* numerical factorization phase */
623329ec9b3SHong Zhang   /*-------------------------------*/
6243d472b54SHong Zhang   lu->id.job = JOB_FACTNUMERIC;
625958c9bccSBarry Smith   if(!lu->id.ICNTL(18)) {
626a7aca84bSHong Zhang     if (!lu->myid) {
627397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
628397b6df1SKris Buschelman       lu->id.a = (mumps_double_complex*)lu->val;
629397b6df1SKris Buschelman #else
630397b6df1SKris Buschelman       lu->id.a = lu->val;
631397b6df1SKris Buschelman #endif
632397b6df1SKris Buschelman     }
633397b6df1SKris Buschelman   } else {
634397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
635397b6df1SKris Buschelman     lu->id.a_loc = (mumps_double_complex*)lu->val;
636397b6df1SKris Buschelman #else
637397b6df1SKris Buschelman     lu->id.a_loc = lu->val;
638397b6df1SKris Buschelman #endif
639397b6df1SKris Buschelman   }
640397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
641397b6df1SKris Buschelman   zmumps_c(&lu->id);
642397b6df1SKris Buschelman #else
643397b6df1SKris Buschelman   dmumps_c(&lu->id);
644397b6df1SKris Buschelman #endif
645397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
64665e19b50SBarry 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));
64765e19b50SBarry 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));
648397b6df1SKris Buschelman   }
64965e19b50SBarry 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));
650397b6df1SKris Buschelman 
6518ada1bb4SHong Zhang   if (lu->size > 1){
65216ebf90aSShri Abhyankar     ierr = PetscTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr);
653a214ac2aSShri Abhyankar     if(isMPIAIJ) {
654dcd589f8SShri Abhyankar       F_diag = ((Mat_MPIAIJ *)(F)->data)->A;
655e09efc27SHong Zhang     } else {
656dcd589f8SShri Abhyankar       F_diag = ((Mat_MPISBAIJ *)(F)->data)->A;
657e09efc27SHong Zhang     }
658e09efc27SHong Zhang     F_diag->assembled = PETSC_TRUE;
659329ec9b3SHong Zhang     if (lu->nSolve){
660329ec9b3SHong Zhang       ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);
6610e83c824SBarry Smith       ierr = PetscFree2(lu->id.sol_loc,lu->id.isol_loc);CHKERRQ(ierr);
662329ec9b3SHong Zhang       ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);
663329ec9b3SHong Zhang     }
6648ada1bb4SHong Zhang   }
665dcd589f8SShri Abhyankar   (F)->assembled   = PETSC_TRUE;
666397b6df1SKris Buschelman   lu->matstruc     = SAME_NONZERO_PATTERN;
667ace87b0dSHong Zhang   lu->CleanUpMUMPS = PETSC_TRUE;
668329ec9b3SHong Zhang   lu->nSolve       = 0;
66967877ebaSShri Abhyankar 
67067877ebaSShri Abhyankar   if (lu->size > 1){
67167877ebaSShri Abhyankar     /* distributed solution */
67267877ebaSShri Abhyankar     lu->id.ICNTL(21) = 1;
67367877ebaSShri Abhyankar     if (!lu->nSolve){
67467877ebaSShri Abhyankar       /* Create x_seq=sol_loc for repeated use */
67567877ebaSShri Abhyankar       PetscInt    lsol_loc;
67667877ebaSShri Abhyankar       PetscScalar *sol_loc;
67767877ebaSShri Abhyankar       lsol_loc = lu->id.INFO(23); /* length of sol_loc */
67867877ebaSShri Abhyankar       ierr = PetscMalloc2(lsol_loc,PetscScalar,&sol_loc,lsol_loc,PetscInt,&lu->id.isol_loc);CHKERRQ(ierr);
67967877ebaSShri Abhyankar       lu->id.lsol_loc = lsol_loc;
68067877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
68167877ebaSShri Abhyankar       lu->id.sol_loc  = (mumps_double_complex*)sol_loc;
68267877ebaSShri Abhyankar #else
68367877ebaSShri Abhyankar       lu->id.sol_loc  = sol_loc;
68467877ebaSShri Abhyankar #endif
68567877ebaSShri Abhyankar       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,lsol_loc,sol_loc,&lu->x_seq);CHKERRQ(ierr);
68667877ebaSShri Abhyankar     }
68767877ebaSShri Abhyankar   }
688397b6df1SKris Buschelman   PetscFunctionReturn(0);
689397b6df1SKris Buschelman }
690397b6df1SKris Buschelman 
691dcd589f8SShri Abhyankar #undef __FUNCT__
692dcd589f8SShri Abhyankar #define __FUNCT__ "PetscSetMUMPSOptions"
693dcd589f8SShri Abhyankar PetscErrorCode PetscSetMUMPSOptions(Mat F, Mat A)
694dcd589f8SShri Abhyankar {
695dcd589f8SShri Abhyankar   Mat_MUMPS        *lu = (Mat_MUMPS*)F->spptr;
696dcd589f8SShri Abhyankar   PetscErrorCode   ierr;
697dcd589f8SShri Abhyankar   PetscInt         icntl;
698dcd589f8SShri Abhyankar   PetscTruth       flg;
699dcd589f8SShri Abhyankar 
700dcd589f8SShri Abhyankar   PetscFunctionBegin;
701dcd589f8SShri Abhyankar   ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
702cb828f0fSHong Zhang   ierr = PetscOptionsTruth("-mat_mumps_view","View MUMPS parameters","None",lu->mumpsview,&lu->mumpsview,PETSC_NULL);CHKERRQ(ierr);
703dcd589f8SShri Abhyankar   if (lu->size == 1){
704dcd589f8SShri Abhyankar     lu->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
705dcd589f8SShri Abhyankar   } else {
706dcd589f8SShri Abhyankar     lu->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
707dcd589f8SShri Abhyankar   }
708dcd589f8SShri Abhyankar 
709dcd589f8SShri Abhyankar   icntl=-1;
710dcd589f8SShri Abhyankar   lu->id.ICNTL(4) = 0;  /* level of printing; overwrite mumps default ICNTL(4)=2 */
711dcd589f8SShri Abhyankar   ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
712dcd589f8SShri Abhyankar   if ((flg && icntl > 0) || PetscLogPrintInfo) {
713dcd589f8SShri Abhyankar     lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */
714dcd589f8SShri Abhyankar   } else { /* no output */
715dcd589f8SShri Abhyankar     lu->id.ICNTL(1) = 0;  /* error message, default= 6 */
716dcd589f8SShri Abhyankar     lu->id.ICNTL(2) = 0;  /* output stream for diagnostic printing, statistics, and warning. default=0 */
717dcd589f8SShri Abhyankar     lu->id.ICNTL(3) = 0; /* output stream for global information, default=6 */
718dcd589f8SShri Abhyankar   }
719dcd589f8SShri 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);
720dcd589f8SShri Abhyankar   icntl=-1;
721dcd589f8SShri Abhyankar   ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
722dcd589f8SShri Abhyankar   if (flg) {
723dcd589f8SShri Abhyankar     if (icntl== 1){
724e32f2f54SBarry 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");
725dcd589f8SShri Abhyankar     } else {
726dcd589f8SShri Abhyankar       lu->id.ICNTL(7) = icntl;
727dcd589f8SShri Abhyankar     }
728dcd589f8SShri Abhyankar   }
729dcd589f8SShri 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);
730dcd589f8SShri 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);
731dcd589f8SShri 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);
732dcd589f8SShri 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);
733dcd589f8SShri 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);
734dcd589f8SShri 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);
735dcd589f8SShri 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);
736dcd589f8SShri Abhyankar   ierr = PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): Schur complement","None",lu->id.ICNTL(19),&lu->id.ICNTL(19),PETSC_NULL);CHKERRQ(ierr);
737dcd589f8SShri Abhyankar 
738dcd589f8SShri 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);
739dcd589f8SShri 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);
740dcd589f8SShri 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);
741dcd589f8SShri 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);
742dcd589f8SShri 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);
743dcd589f8SShri Abhyankar   ierr = PetscOptionsInt("-mat_mumps_icntl_27","ICNTL(27): experimental parameter","None",lu->id.ICNTL(27),&lu->id.ICNTL(27),PETSC_NULL);CHKERRQ(ierr);
744dcd589f8SShri Abhyankar 
745dcd589f8SShri 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);
746dcd589f8SShri 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);
747dcd589f8SShri 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);
748dcd589f8SShri 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);
749dcd589f8SShri 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);
750dcd589f8SShri Abhyankar   PetscOptionsEnd();
751dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
752dcd589f8SShri Abhyankar }
753dcd589f8SShri Abhyankar 
754dcd589f8SShri Abhyankar #undef __FUNCT__
755dcd589f8SShri Abhyankar #define __FUNCT__ "PetscInitializeMUMPS"
756f697e70eSHong Zhang PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS* mumps)
757dcd589f8SShri Abhyankar {
758dcd589f8SShri Abhyankar   PetscErrorCode  ierr;
759dcd589f8SShri Abhyankar 
760dcd589f8SShri Abhyankar   PetscFunctionBegin;
761f697e70eSHong Zhang   ierr = MPI_Comm_rank(((PetscObject)A)->comm, &mumps->myid);
762f697e70eSHong Zhang   ierr = MPI_Comm_size(((PetscObject)A)->comm,&mumps->size);CHKERRQ(ierr);
763f697e70eSHong Zhang   ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(mumps->comm_mumps));CHKERRQ(ierr);
764f697e70eSHong Zhang   mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps);
765f697e70eSHong Zhang 
766f697e70eSHong Zhang   mumps->id.job = JOB_INIT;
767f697e70eSHong Zhang   mumps->id.par = 1;  /* host participates factorizaton and solve */
768f697e70eSHong Zhang   mumps->id.sym = mumps->sym;
769dcd589f8SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
770f697e70eSHong Zhang   zmumps_c(&mumps->id);
771dcd589f8SShri Abhyankar #else
772f697e70eSHong Zhang   dmumps_c(&mumps->id);
773dcd589f8SShri Abhyankar #endif
774f697e70eSHong Zhang 
775f697e70eSHong Zhang   mumps->CleanUpMUMPS = PETSC_FALSE;
776f697e70eSHong Zhang   mumps->scat_rhs     = PETSC_NULL;
777f697e70eSHong Zhang   mumps->scat_sol     = PETSC_NULL;
778f697e70eSHong Zhang   mumps->nSolve       = 0;
779dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
780dcd589f8SShri Abhyankar }
781dcd589f8SShri Abhyankar 
782397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */
783397b6df1SKris Buschelman #undef __FUNCT__
784f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
7850481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
786b24902e0SBarry Smith {
787719d5645SBarry Smith   Mat_MUMPS          *lu = (Mat_MUMPS*)F->spptr;
788dcd589f8SShri Abhyankar   PetscErrorCode     ierr;
789bccb9932SShri Abhyankar   MatReuse           reuse;
79067877ebaSShri Abhyankar   Vec                b;
79167877ebaSShri Abhyankar   IS                 is_iden;
79267877ebaSShri Abhyankar   const PetscInt     M = A->rmap->N;
793397b6df1SKris Buschelman 
794397b6df1SKris Buschelman   PetscFunctionBegin;
795b24902e0SBarry Smith   lu->matstruc = DIFFERENT_NONZERO_PATTERN;
796dcd589f8SShri Abhyankar 
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   /*----------------*/
8053d472b54SHong Zhang   lu->id.job = JOB_FACTSYMBOLIC;
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->matstruc = DIFFERENT_NONZERO_PATTERN;
874dcd589f8SShri Abhyankar 
875dcd589f8SShri Abhyankar   /* Set MUMPS options */
876dcd589f8SShri Abhyankar   ierr = PetscSetMUMPSOptions(F,A);CHKERRQ(ierr);
877dcd589f8SShri Abhyankar 
878bccb9932SShri Abhyankar   reuse = MAT_INITIAL_MATRIX;
879bccb9932SShri Abhyankar   ierr = (*lu->ConvertToTriples)(A, 1, reuse, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
88067877ebaSShri Abhyankar 
88167877ebaSShri Abhyankar   /* analysis phase */
88267877ebaSShri Abhyankar   /*----------------*/
8833d472b54SHong Zhang   lu->id.job = JOB_FACTSYMBOLIC;
88467877ebaSShri Abhyankar   lu->id.n = M;
88567877ebaSShri Abhyankar   switch (lu->id.ICNTL(18)){
88667877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
88767877ebaSShri Abhyankar     if (!lu->myid) {
88867877ebaSShri Abhyankar       lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
88967877ebaSShri Abhyankar       if (lu->id.ICNTL(6)>1){
89067877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
89167877ebaSShri Abhyankar         lu->id.a = (mumps_double_complex*)lu->val;
89267877ebaSShri Abhyankar #else
89367877ebaSShri Abhyankar         lu->id.a = lu->val;
89467877ebaSShri Abhyankar #endif
89567877ebaSShri Abhyankar       }
89667877ebaSShri Abhyankar     }
89767877ebaSShri Abhyankar     break;
89867877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
89967877ebaSShri Abhyankar     lu->id.nz_loc = lu->nz;
90067877ebaSShri Abhyankar     lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
90167877ebaSShri Abhyankar     if (lu->id.ICNTL(6)>1) {
90267877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
90367877ebaSShri Abhyankar       lu->id.a_loc = (mumps_double_complex*)lu->val;
90467877ebaSShri Abhyankar #else
90567877ebaSShri Abhyankar       lu->id.a_loc = lu->val;
90667877ebaSShri Abhyankar #endif
90767877ebaSShri Abhyankar     }
90867877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
90967877ebaSShri Abhyankar     if (!lu->myid){
91067877ebaSShri Abhyankar       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr);
91167877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
91267877ebaSShri Abhyankar     } else {
91367877ebaSShri Abhyankar       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr);
91467877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
91567877ebaSShri Abhyankar     }
91667877ebaSShri Abhyankar     ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr);
91767877ebaSShri Abhyankar     ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
91867877ebaSShri Abhyankar     ierr = VecSetFromOptions(b);CHKERRQ(ierr);
91967877ebaSShri Abhyankar 
92067877ebaSShri Abhyankar     ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr);
92167877ebaSShri Abhyankar     ierr = ISDestroy(is_iden);CHKERRQ(ierr);
92267877ebaSShri Abhyankar     ierr = VecDestroy(b);CHKERRQ(ierr);
92367877ebaSShri Abhyankar     break;
92467877ebaSShri Abhyankar     }
92567877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
92667877ebaSShri Abhyankar   zmumps_c(&lu->id);
92767877ebaSShri Abhyankar #else
92867877ebaSShri Abhyankar   dmumps_c(&lu->id);
92967877ebaSShri Abhyankar #endif
93067877ebaSShri 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));
93167877ebaSShri Abhyankar 
932450b117fSShri Abhyankar   F->ops->lufactornumeric  = MatFactorNumeric_MUMPS;
933dcd589f8SShri Abhyankar   F->ops->solve            = MatSolve_MUMPS;
934450b117fSShri Abhyankar   PetscFunctionReturn(0);
935450b117fSShri Abhyankar }
936b24902e0SBarry Smith 
937397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */
938397b6df1SKris Buschelman #undef __FUNCT__
93967877ebaSShri Abhyankar #define __FUNCT__ "MatCholeskyFactorSymbolic_MUMPS"
94067877ebaSShri Abhyankar PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
941b24902e0SBarry Smith {
9422792810eSHong Zhang   Mat_MUMPS          *lu = (Mat_MUMPS*)F->spptr;
943dcd589f8SShri Abhyankar   PetscErrorCode     ierr;
944bccb9932SShri Abhyankar   MatReuse           reuse;
94567877ebaSShri Abhyankar   Vec                b;
94667877ebaSShri Abhyankar   IS                 is_iden;
94767877ebaSShri Abhyankar   const PetscInt     M = A->rmap->N;
948397b6df1SKris Buschelman 
949397b6df1SKris Buschelman   PetscFunctionBegin;
950b24902e0SBarry Smith   lu->matstruc = DIFFERENT_NONZERO_PATTERN;
951dcd589f8SShri Abhyankar 
952dcd589f8SShri Abhyankar   /* Set MUMPS options */
953dcd589f8SShri Abhyankar   ierr = PetscSetMUMPSOptions(F,A);CHKERRQ(ierr);
954dcd589f8SShri Abhyankar 
955bccb9932SShri Abhyankar   reuse = MAT_INITIAL_MATRIX;
956bccb9932SShri Abhyankar   ierr = (*lu->ConvertToTriples)(A, 1 , reuse, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
957dcd589f8SShri Abhyankar 
95867877ebaSShri Abhyankar   /* analysis phase */
95967877ebaSShri Abhyankar   /*----------------*/
9603d472b54SHong Zhang   lu->id.job = JOB_FACTSYMBOLIC;
96167877ebaSShri Abhyankar   lu->id.n = M;
96267877ebaSShri Abhyankar   switch (lu->id.ICNTL(18)){
96367877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
96467877ebaSShri Abhyankar     if (!lu->myid) {
96567877ebaSShri Abhyankar       lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
96667877ebaSShri Abhyankar       if (lu->id.ICNTL(6)>1){
96767877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
96867877ebaSShri Abhyankar         lu->id.a = (mumps_double_complex*)lu->val;
96967877ebaSShri Abhyankar #else
97067877ebaSShri Abhyankar         lu->id.a = lu->val;
97167877ebaSShri Abhyankar #endif
97267877ebaSShri Abhyankar       }
97367877ebaSShri Abhyankar     }
97467877ebaSShri Abhyankar     break;
97567877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
97667877ebaSShri Abhyankar     lu->id.nz_loc = lu->nz;
97767877ebaSShri Abhyankar     lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
97867877ebaSShri Abhyankar     if (lu->id.ICNTL(6)>1) {
97967877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
98067877ebaSShri Abhyankar       lu->id.a_loc = (mumps_double_complex*)lu->val;
98167877ebaSShri Abhyankar #else
98267877ebaSShri Abhyankar       lu->id.a_loc = lu->val;
98367877ebaSShri Abhyankar #endif
98467877ebaSShri Abhyankar     }
98567877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
98667877ebaSShri Abhyankar     if (!lu->myid){
98767877ebaSShri Abhyankar       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr);
98867877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
98967877ebaSShri Abhyankar     } else {
99067877ebaSShri Abhyankar       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr);
99167877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
99267877ebaSShri Abhyankar     }
99367877ebaSShri Abhyankar     ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr);
99467877ebaSShri Abhyankar     ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
99567877ebaSShri Abhyankar     ierr = VecSetFromOptions(b);CHKERRQ(ierr);
99667877ebaSShri Abhyankar 
99767877ebaSShri Abhyankar     ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr);
99867877ebaSShri Abhyankar     ierr = ISDestroy(is_iden);CHKERRQ(ierr);
99967877ebaSShri Abhyankar     ierr = VecDestroy(b);CHKERRQ(ierr);
100067877ebaSShri Abhyankar     break;
100167877ebaSShri Abhyankar     }
100267877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
100367877ebaSShri Abhyankar   zmumps_c(&lu->id);
100467877ebaSShri Abhyankar #else
100567877ebaSShri Abhyankar   dmumps_c(&lu->id);
100667877ebaSShri Abhyankar #endif
100767877ebaSShri 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));
100867877ebaSShri Abhyankar 
10092792810eSHong Zhang   F->ops->choleskyfactornumeric =  MatFactorNumeric_MUMPS;
1010dcd589f8SShri Abhyankar   F->ops->solve                 =  MatSolve_MUMPS;
1011db4efbfdSBarry Smith #if !defined(PETSC_USE_COMPLEX)
1012dcd589f8SShri Abhyankar   (F)->ops->getinertia          =  MatGetInertia_SBAIJMUMPS;
1013db4efbfdSBarry Smith #endif
1014b24902e0SBarry Smith   PetscFunctionReturn(0);
1015b24902e0SBarry Smith }
1016b24902e0SBarry Smith 
1017397b6df1SKris Buschelman #undef __FUNCT__
1018f6c57405SHong Zhang #define __FUNCT__ "MatFactorInfo_MUMPS"
101974ed9c26SBarry Smith PetscErrorCode MatFactorInfo_MUMPS(Mat A,PetscViewer viewer)
102074ed9c26SBarry Smith {
1021f6c57405SHong Zhang   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
1022f6c57405SHong Zhang   PetscErrorCode ierr;
1023f6c57405SHong Zhang 
1024f6c57405SHong Zhang   PetscFunctionBegin;
1025f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):        %d \n",lu->id.ICNTL(1));CHKERRQ(ierr);
1026f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg):%d \n",lu->id.ICNTL(2));CHKERRQ(ierr);
1027f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):  %d \n",lu->id.ICNTL(3));CHKERRQ(ierr);
1028f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):       %d \n",lu->id.ICNTL(4));CHKERRQ(ierr);
1029f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):        %d \n",lu->id.ICNTL(5));CHKERRQ(ierr);
1030f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):       %d \n",lu->id.ICNTL(6));CHKERRQ(ierr);
1031f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (matrix ordering):         %d \n",lu->id.ICNTL(7));CHKERRQ(ierr);
1032f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scalling strategy):       %d \n",lu->id.ICNTL(8));CHKERRQ(ierr);
1033f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(9) (A/A^T x=b is solved):     %d \n",lu->id.ICNTL(9));CHKERRQ(ierr);
1034f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr);
1035f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):         %d \n",lu->id.ICNTL(11));CHKERRQ(ierr);
103634ed7027SBarry Smith   if (lu->id.ICNTL(11)>0) {
103734ed7027SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"        RINFOG(4) (inf norm of input mat):        %g\n",lu->id.RINFOG(4));CHKERRQ(ierr);
103834ed7027SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"        RINFOG(5) (inf norm of solution):         %g\n",lu->id.RINFOG(5));CHKERRQ(ierr);
103934ed7027SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"        RINFOG(6) (inf norm of residual):         %g\n",lu->id.RINFOG(6));CHKERRQ(ierr);
104034ed7027SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"        RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",lu->id.RINFOG(7),lu->id.RINFOG(8));CHKERRQ(ierr);
104134ed7027SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"        RINFOG(9) (error estimate):               %g \n",lu->id.RINFOG(9));CHKERRQ(ierr);
104234ed7027SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"        RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",lu->id.RINFOG(10),lu->id.RINFOG(11));CHKERRQ(ierr);
1043f6c57405SHong Zhang 
1044f6c57405SHong Zhang   }
1045f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",lu->id.ICNTL(12));CHKERRQ(ierr);
1046f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",lu->id.ICNTL(13));CHKERRQ(ierr);
1047f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr);
1048f6c57405SHong Zhang   /* ICNTL(15-17) not used */
1049f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",lu->id.ICNTL(18));CHKERRQ(ierr);
1050f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Shur complement info):                       %d \n",lu->id.ICNTL(19));CHKERRQ(ierr);
1051f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",lu->id.ICNTL(20));CHKERRQ(ierr);
1052f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d \n",lu->id.ICNTL(21));CHKERRQ(ierr);
1053c0165424SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",lu->id.ICNTL(22));CHKERRQ(ierr);
1054c0165424SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",lu->id.ICNTL(23));CHKERRQ(ierr);
1055c0165424SHong Zhang 
1056c0165424SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",lu->id.ICNTL(24));CHKERRQ(ierr);
1057c0165424SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",lu->id.ICNTL(25));CHKERRQ(ierr);
1058c0165424SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",lu->id.ICNTL(26));CHKERRQ(ierr);
1059c0165424SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",lu->id.ICNTL(27));CHKERRQ(ierr);
1060f6c57405SHong Zhang 
1061f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",lu->id.CNTL(1));CHKERRQ(ierr);
1062f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr);
1063f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",lu->id.CNTL(3));CHKERRQ(ierr);
1064f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",lu->id.CNTL(4));CHKERRQ(ierr);
1065c0165424SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",lu->id.CNTL(5));CHKERRQ(ierr);
1066f6c57405SHong Zhang 
1067f6c57405SHong Zhang   /* infomation local to each processor */
106834ed7027SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer, "              RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);
106934ed7027SBarry Smith   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"              [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr);
107034ed7027SBarry Smith   ierr = PetscViewerFlush(viewer);
107134ed7027SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer, "              RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);
107234ed7027SBarry Smith   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"              [%d]  %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr);
107334ed7027SBarry Smith   ierr = PetscViewerFlush(viewer);
107434ed7027SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer, "              RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);
107534ed7027SBarry Smith   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"              [%d]  %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr);
107634ed7027SBarry Smith   ierr = PetscViewerFlush(viewer);
1077f6c57405SHong Zhang 
107834ed7027SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer, "              INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr);
107934ed7027SBarry Smith   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"              [%d] %d \n",lu->myid,lu->id.INFO(15));CHKERRQ(ierr);
108034ed7027SBarry Smith   ierr = PetscViewerFlush(viewer);
1081f6c57405SHong Zhang 
108234ed7027SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer, "              INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr);
108334ed7027SBarry Smith   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"              [%d] %d \n",lu->myid,lu->id.INFO(16));CHKERRQ(ierr);
108434ed7027SBarry Smith   ierr = PetscViewerFlush(viewer);
1085f6c57405SHong Zhang 
108634ed7027SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer, "              INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);
108734ed7027SBarry Smith   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"              [%d] %d \n",lu->myid,lu->id.INFO(23));CHKERRQ(ierr);
108834ed7027SBarry Smith   ierr = PetscViewerFlush(viewer);
1089f6c57405SHong Zhang 
1090f6c57405SHong Zhang   if (!lu->myid){ /* information from the host */
1091f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr);
1092f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr);
1093f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr);
1094f6c57405SHong Zhang 
1095f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr);
1096f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr);
1097f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr);
1098f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr);
1099f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr);
1100f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr);
1101f6c57405SHong 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);
1102f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",lu->id.INFOG(10));CHKERRQ(ierr);
1103f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",lu->id.INFOG(11));CHKERRQ(ierr);
1104f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr);
1105f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr);
1106f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr);
1107f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr);
1108f6c57405SHong 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);
1109f6c57405SHong 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);
1110f6c57405SHong 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);
1111f6c57405SHong 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);
1112f6c57405SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr);
1113f6c57405SHong 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);
1114f6c57405SHong 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);
1115f6c57405SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",lu->id.INFOG(23));CHKERRQ(ierr);
1116f6c57405SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",lu->id.INFOG(24));CHKERRQ(ierr);
1117f6c57405SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",lu->id.INFOG(25));CHKERRQ(ierr);
1118f6c57405SHong Zhang   }
1119f6c57405SHong Zhang   PetscFunctionReturn(0);
1120f6c57405SHong Zhang }
1121f6c57405SHong Zhang 
1122f6c57405SHong Zhang #undef __FUNCT__
1123f6c57405SHong Zhang #define __FUNCT__ "MatView_MUMPS"
1124b24902e0SBarry Smith PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
1125b24902e0SBarry Smith {
1126f6c57405SHong Zhang   PetscErrorCode    ierr;
1127f6c57405SHong Zhang   PetscTruth        iascii;
1128f6c57405SHong Zhang   PetscViewerFormat format;
1129cb828f0fSHong Zhang   Mat_MUMPS         *mumps=(Mat_MUMPS*)A->spptr;
1130f6c57405SHong Zhang 
1131f6c57405SHong Zhang   PetscFunctionBegin;
1132cb828f0fSHong Zhang   /* check if matrix is mumps type */
1133cb828f0fSHong Zhang   if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0);
1134cb828f0fSHong Zhang 
11352692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1136f6c57405SHong Zhang   if (iascii) {
1137f6c57405SHong Zhang     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1138cb828f0fSHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO || mumps->mumpsview){
1139cb828f0fSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
1140cb828f0fSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                  %d \n",mumps->id.sym);CHKERRQ(ierr);
1141cb828f0fSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):           %d \n",mumps->id.par);CHKERRQ(ierr);
1142cb828f0fSHong Zhang       if (mumps->mumpsview){ /* View all MUMPS parameters */
1143f6c57405SHong Zhang         ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr);
1144f6c57405SHong Zhang       }
1145f6c57405SHong Zhang     }
1146cb828f0fSHong Zhang   }
1147f6c57405SHong Zhang   PetscFunctionReturn(0);
1148f6c57405SHong Zhang }
1149f6c57405SHong Zhang 
115035bd34faSBarry Smith #undef __FUNCT__
115135bd34faSBarry Smith #define __FUNCT__ "MatGetInfo_MUMPS"
115235bd34faSBarry Smith PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
115335bd34faSBarry Smith {
1154cb828f0fSHong Zhang   Mat_MUMPS  *mumps =(Mat_MUMPS*)A->spptr;
115535bd34faSBarry Smith 
115635bd34faSBarry Smith   PetscFunctionBegin;
115735bd34faSBarry Smith   info->block_size        = 1.0;
1158cb828f0fSHong Zhang   info->nz_allocated      = mumps->id.INFOG(20);
1159cb828f0fSHong Zhang   info->nz_used           = mumps->id.INFOG(20);
116035bd34faSBarry Smith   info->nz_unneeded       = 0.0;
116135bd34faSBarry Smith   info->assemblies        = 0.0;
116235bd34faSBarry Smith   info->mallocs           = 0.0;
116335bd34faSBarry Smith   info->memory            = 0.0;
116435bd34faSBarry Smith   info->fill_ratio_given  = 0;
116535bd34faSBarry Smith   info->fill_ratio_needed = 0;
116635bd34faSBarry Smith   info->factor_mallocs    = 0;
116735bd34faSBarry Smith   PetscFunctionReturn(0);
116835bd34faSBarry Smith }
116935bd34faSBarry Smith 
117024b6179bSKris Buschelman /*MC
11712692d6eeSBarry Smith   MATSOLVERMUMPS -  A matrix type providing direct solvers (LU and Cholesky) for
117224b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
117324b6179bSKris Buschelman 
117441c8de11SBarry Smith   Works with MATAIJ and MATSBAIJ matrices
117524b6179bSKris Buschelman 
117624b6179bSKris Buschelman   Options Database Keys:
1177*fb8376fbSHong Zhang + -mat_mumps_icntl_4 <0,...,4> - print level
117824b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
117924b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
118024b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
118124b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
118294b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
118324b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
118424b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
118524b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
118624b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
118724b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
118824b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
118924b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
119024b6179bSKris Buschelman 
119124b6179bSKris Buschelman   Level: beginner
119224b6179bSKris Buschelman 
119341c8de11SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage
119441c8de11SBarry Smith 
119524b6179bSKris Buschelman M*/
119624b6179bSKris Buschelman 
11972877fffaSHong Zhang EXTERN_C_BEGIN
119835bd34faSBarry Smith #undef __FUNCT__
119935bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_mumps"
120035bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type)
120135bd34faSBarry Smith {
120235bd34faSBarry Smith   PetscFunctionBegin;
12032692d6eeSBarry Smith   *type = MATSOLVERMUMPS;
120435bd34faSBarry Smith   PetscFunctionReturn(0);
120535bd34faSBarry Smith }
120635bd34faSBarry Smith EXTERN_C_END
120735bd34faSBarry Smith 
120835bd34faSBarry Smith EXTERN_C_BEGIN
1209bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI AIJ matrices */
12102877fffaSHong Zhang #undef __FUNCT__
1211bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_aij_mumps"
1212bccb9932SShri Abhyankar PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
12132877fffaSHong Zhang {
12142877fffaSHong Zhang   Mat            B;
12152877fffaSHong Zhang   PetscErrorCode ierr;
12162877fffaSHong Zhang   Mat_MUMPS      *mumps;
1217bccb9932SShri Abhyankar   PetscTruth     isSeqAIJ;
12182877fffaSHong Zhang 
12192877fffaSHong Zhang   PetscFunctionBegin;
12202877fffaSHong Zhang   /* Create the factorization matrix */
1221bccb9932SShri Abhyankar   ierr = PetscTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr);
12222877fffaSHong Zhang   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
12232877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
12242877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1225bccb9932SShri Abhyankar   if (isSeqAIJ) {
12262877fffaSHong Zhang     ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
1227bccb9932SShri Abhyankar   } else {
1228bccb9932SShri Abhyankar     ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
1229bccb9932SShri Abhyankar   }
12302877fffaSHong Zhang 
123116ebf90aSShri Abhyankar   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
12322877fffaSHong Zhang   B->ops->view             = MatView_MUMPS;
123335bd34faSBarry Smith   B->ops->getinfo          = MatGetInfo_MUMPS;
123435bd34faSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
1235a1f19f5aSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSetMumpsIcntl_C","MatSetMumpsIcntl",MatSetMumpsIcntl);CHKERRQ(ierr);
1236450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
1237450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
1238d5f3da31SBarry Smith     B->factortype = MAT_FACTOR_LU;
1239bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
1240bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
1241746480a1SHong Zhang     mumps->sym = 0;
1242dcd589f8SShri Abhyankar   } else {
124367877ebaSShri Abhyankar     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
1244450b117fSShri Abhyankar     B->factortype = MAT_FACTOR_CHOLESKY;
1245bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
1246bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
12476fdc2a6dSBarry Smith     if (A->spd_set && A->spd) mumps->sym = 1;
12486fdc2a6dSBarry Smith     else                      mumps->sym = 2;
1249450b117fSShri Abhyankar   }
12502877fffaSHong Zhang 
12512877fffaSHong Zhang   mumps->isAIJ        = PETSC_TRUE;
12522877fffaSHong Zhang   mumps->MatDestroy   = B->ops->destroy;
12532877fffaSHong Zhang   B->ops->destroy     = MatDestroy_MUMPS;
12542877fffaSHong Zhang   B->spptr            = (void*)mumps;
1255f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1256746480a1SHong Zhang 
12572877fffaSHong Zhang   *F = B;
12582877fffaSHong Zhang   PetscFunctionReturn(0);
12592877fffaSHong Zhang }
12602877fffaSHong Zhang EXTERN_C_END
12612877fffaSHong Zhang 
1262bccb9932SShri Abhyankar 
12632877fffaSHong Zhang EXTERN_C_BEGIN
1264bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI SBAIJ matrices */
12652877fffaSHong Zhang #undef __FUNCT__
1266bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_sbaij_mumps"
1267bccb9932SShri Abhyankar PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
12682877fffaSHong Zhang {
12692877fffaSHong Zhang   Mat            B;
12702877fffaSHong Zhang   PetscErrorCode ierr;
12712877fffaSHong Zhang   Mat_MUMPS      *mumps;
1272bccb9932SShri Abhyankar   PetscTruth     isSeqSBAIJ;
12732877fffaSHong Zhang 
12742877fffaSHong Zhang   PetscFunctionBegin;
1275e7e72b3dSBarry Smith   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
1276bccb9932SShri 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");
1277bccb9932SShri Abhyankar   ierr = PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
12782877fffaSHong Zhang   /* Create the factorization matrix */
12792877fffaSHong Zhang   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
12802877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
12812877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
128216ebf90aSShri Abhyankar   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
1283bccb9932SShri Abhyankar   if (isSeqSBAIJ) {
1284bccb9932SShri Abhyankar     ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
128516ebf90aSShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
1286dcd589f8SShri Abhyankar   } else {
1287bccb9932SShri Abhyankar     ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
1288bccb9932SShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
1289bccb9932SShri Abhyankar   }
1290bccb9932SShri Abhyankar 
129167877ebaSShri Abhyankar   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
1292bccb9932SShri Abhyankar   B->ops->view                   = MatView_MUMPS;
1293bccb9932SShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
1294a1f19f5aSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSetMumpsIcntl_C","MatSetMumpsIcntl",MatSetMumpsIcntl);CHKERRQ(ierr);
1295f4762488SHong Zhang   B->factortype                  = MAT_FACTOR_CHOLESKY;
12966fdc2a6dSBarry Smith   if (A->spd_set && A->spd) mumps->sym = 1;
12976fdc2a6dSBarry Smith   else                      mumps->sym = 2;
1298a214ac2aSShri Abhyankar 
1299bccb9932SShri Abhyankar   mumps->isAIJ        = PETSC_FALSE;
1300f3c0ef26SHong Zhang   mumps->MatDestroy   = B->ops->destroy;
1301f3c0ef26SHong Zhang   B->ops->destroy     = MatDestroy_MUMPS;
13022877fffaSHong Zhang   B->spptr            = (void*)mumps;
1303f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1304746480a1SHong Zhang 
13052877fffaSHong Zhang   *F = B;
13062877fffaSHong Zhang   PetscFunctionReturn(0);
13072877fffaSHong Zhang }
13082877fffaSHong Zhang EXTERN_C_END
130997969023SHong Zhang 
1310450b117fSShri Abhyankar EXTERN_C_BEGIN
1311450b117fSShri Abhyankar #undef __FUNCT__
1312bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_baij_mumps"
1313bccb9932SShri Abhyankar PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
131467877ebaSShri Abhyankar {
131567877ebaSShri Abhyankar   Mat            B;
131667877ebaSShri Abhyankar   PetscErrorCode ierr;
131767877ebaSShri Abhyankar   Mat_MUMPS      *mumps;
1318bccb9932SShri Abhyankar   PetscTruth     isSeqBAIJ;
131967877ebaSShri Abhyankar 
132067877ebaSShri Abhyankar   PetscFunctionBegin;
132167877ebaSShri Abhyankar   /* Create the factorization matrix */
1322bccb9932SShri Abhyankar   ierr = PetscTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr);
132367877ebaSShri Abhyankar   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
132467877ebaSShri Abhyankar   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
132567877ebaSShri Abhyankar   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1326bccb9932SShri Abhyankar   if (isSeqBAIJ) {
132767877ebaSShri Abhyankar     ierr = MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,PETSC_NULL);CHKERRQ(ierr);
1328bccb9932SShri Abhyankar   } else {
132967877ebaSShri Abhyankar     ierr = MatMPIBAIJSetPreallocation(B,A->rmap->bs,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
1330bccb9932SShri Abhyankar   }
1331450b117fSShri Abhyankar 
133267877ebaSShri Abhyankar   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
1333450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
1334450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
1335450b117fSShri Abhyankar     B->factortype = MAT_FACTOR_LU;
1336bccb9932SShri Abhyankar     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
1337bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
1338746480a1SHong Zhang     mumps->sym = 0;
1339746480a1SHong Zhang   } else {
1340746480a1SHong Zhang     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n");
1341450b117fSShri Abhyankar   }
1342bccb9932SShri Abhyankar 
1343450b117fSShri Abhyankar   B->ops->view             = MatView_MUMPS;
1344450b117fSShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
1345a1f19f5aSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSetMumpsIcntl_C","MatSetMumpsIcntl",MatSetMumpsIcntl);CHKERRQ(ierr);
1346450b117fSShri Abhyankar 
1347450b117fSShri Abhyankar   mumps->isAIJ        = PETSC_TRUE;
1348450b117fSShri Abhyankar   mumps->MatDestroy   = B->ops->destroy;
1349450b117fSShri Abhyankar   B->ops->destroy     = MatDestroy_MUMPS;
1350450b117fSShri Abhyankar   B->spptr            = (void*)mumps;
1351f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1352746480a1SHong Zhang 
1353450b117fSShri Abhyankar   *F = B;
1354450b117fSShri Abhyankar   PetscFunctionReturn(0);
1355450b117fSShri Abhyankar }
1356450b117fSShri Abhyankar EXTERN_C_END
1357a214ac2aSShri Abhyankar 
135861288eaaSHong Zhang /* -------------------------------------------------------------------------------------------*/
1359ed4114aaSSatish Balay #undef __FUNCT__
13608b1197c0SBarry Smith #define __FUNCT__ "MatSetMumpsIcntl"
136161288eaaSHong Zhang /*@
1362a1f19f5aSHong Zhang   MatSetMumpsIcntl - Set MUMPS parameter ICNTL()
136361288eaaSHong Zhang 
13643f9fe445SBarry Smith    Logically Collective on Mat
136561288eaaSHong Zhang 
136661288eaaSHong Zhang    Input Parameters:
136761288eaaSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1368a1f19f5aSHong Zhang .  icntl - index of MUMPS parameter array ICNTL()
1369a1f19f5aSHong Zhang -  ival - value of MUMPS ICNTL(icntl)
137061288eaaSHong Zhang 
137161288eaaSHong Zhang   Options Database:
1372a1f19f5aSHong Zhang .   -mat_mumps_icntl_<icntl> <ival>
137361288eaaSHong Zhang 
137461288eaaSHong Zhang    Level: beginner
137561288eaaSHong Zhang 
137661288eaaSHong Zhang    References: MUMPS Users' Guide
137761288eaaSHong Zhang 
137861288eaaSHong Zhang .seealso: MatGetFactor()
137961288eaaSHong Zhang @*/
1380a1f19f5aSHong Zhang PetscErrorCode MatSetMumpsIcntl(Mat F,PetscInt icntl,PetscInt ival)
138197969023SHong Zhang {
1382dcd589f8SShri Abhyankar   Mat_MUMPS      *lu =(Mat_MUMPS*)(F)->spptr;
138397969023SHong Zhang 
138497969023SHong Zhang   PetscFunctionBegin;
13854068e74aSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
13864068e74aSHong Zhang   PetscValidLogicalCollectiveInt(F,ival,3);
1387a1f19f5aSHong Zhang   lu->id.ICNTL(icntl) = ival;
138897969023SHong Zhang   PetscFunctionReturn(0);
138997969023SHong Zhang }
139097969023SHong Zhang 
1391