xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision 329ec9b3932dd2adb9721d1bf6d0bfd162ea8cbf)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
21c2a3de1SBarry Smith 
3397b6df1SKris Buschelman /*
4c2b5dc30SHong Zhang     Provides an interface to the MUMPS sparse solver
5397b6df1SKris Buschelman */
6397b6df1SKris Buschelman #include "src/mat/impls/aij/seq/aij.h"
7397b6df1SKris Buschelman #include "src/mat/impls/aij/mpi/mpiaij.h"
8397b6df1SKris Buschelman #include "src/mat/impls/sbaij/seq/sbaij.h"
9397b6df1SKris Buschelman #include "src/mat/impls/sbaij/mpi/mpisbaij.h"
10397b6df1SKris Buschelman 
11397b6df1SKris Buschelman EXTERN_C_BEGIN
12397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
13397b6df1SKris Buschelman #include "zmumps_c.h"
14397b6df1SKris Buschelman #else
15397b6df1SKris Buschelman #include "dmumps_c.h"
16397b6df1SKris Buschelman #endif
17397b6df1SKris Buschelman EXTERN_C_END
18397b6df1SKris Buschelman #define JOB_INIT -1
19397b6df1SKris Buschelman #define JOB_END -2
20397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */
21397b6df1SKris Buschelman #define ICNTL(I) icntl[(I)-1]
22397b6df1SKris Buschelman #define CNTL(I) cntl[(I)-1]
23397b6df1SKris Buschelman #define INFOG(I) infog[(I)-1]
24a7aca84bSHong Zhang #define INFO(I) info[(I)-1]
25397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1]
26adc1d99fSHong Zhang #define RINFO(I) rinfo[(I)-1]
27397b6df1SKris Buschelman 
28397b6df1SKris Buschelman typedef struct {
29397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
30397b6df1SKris Buschelman   ZMUMPS_STRUC_C id;
31397b6df1SKris Buschelman #else
32397b6df1SKris Buschelman   DMUMPS_STRUC_C id;
33397b6df1SKris Buschelman #endif
34397b6df1SKris Buschelman   MatStructure   matstruc;
35c1490034SHong Zhang   PetscMPIInt    myid,size;
36*329ec9b3SHong Zhang   PetscInt       *irn,*jcn,sym,nSolve;
37397b6df1SKris Buschelman   PetscScalar    *val;
38397b6df1SKris Buschelman   MPI_Comm       comm_mumps;
39*329ec9b3SHong Zhang   VecScatter     scat_rhs, scat_sol;
40c338a77dSKris Buschelman   PetscTruth     isAIJ,CleanUpMUMPS;
41*329ec9b3SHong Zhang   Vec            b_seq,x_seq;
426849ba73SBarry Smith   PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
436849ba73SBarry Smith   PetscErrorCode (*MatView)(Mat,PetscViewer);
446849ba73SBarry Smith   PetscErrorCode (*MatAssemblyEnd)(Mat,MatAssemblyType);
456849ba73SBarry Smith   PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
466849ba73SBarry Smith   PetscErrorCode (*MatCholeskyFactorSymbolic)(Mat,IS,MatFactorInfo*,Mat*);
476849ba73SBarry Smith   PetscErrorCode (*MatDestroy)(Mat);
486849ba73SBarry Smith   PetscErrorCode (*specialdestroy)(Mat);
496849ba73SBarry Smith   PetscErrorCode (*MatPreallocate)(Mat,int,int,int*,int,int*);
50f0c56d0fSKris Buschelman } Mat_MUMPS;
51f0c56d0fSKris Buschelman 
52dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);
53892f6c3fSKris Buschelman EXTERN_C_BEGIN
5475179d2cSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SBAIJ_SBAIJMUMPS(Mat,MatType,MatReuse,Mat*);
55892f6c3fSKris Buschelman EXTERN_C_END
56397b6df1SKris Buschelman /* convert Petsc mpiaij matrix to triples: row[nz], col[nz], val[nz] */
57397b6df1SKris Buschelman /*
58397b6df1SKris Buschelman   input:
5975747be1SHong Zhang     A       - matrix in mpiaij or mpisbaij (bs=1) format
60397b6df1SKris Buschelman     shift   - 0: C style output triple; 1: Fortran style output triple.
61397b6df1SKris Buschelman     valOnly - FALSE: spaces are allocated and values are set for the triple
62397b6df1SKris Buschelman               TRUE:  only the values in v array are updated
63397b6df1SKris Buschelman   output:
64397b6df1SKris Buschelman     nnz     - dim of r, c, and v (number of local nonzero entries of A)
65397b6df1SKris Buschelman     r, c, v - row and col index, matrix values (matrix triples)
66397b6df1SKris Buschelman  */
67dfbe8321SBarry Smith PetscErrorCode MatConvertToTriples(Mat A,int shift,PetscTruth valOnly,int *nnz,int **r, int **c, PetscScalar **v) {
68c1490034SHong Zhang   PetscInt       *ai, *aj, *bi, *bj, rstart,nz, *garray;
69dfbe8321SBarry Smith   PetscErrorCode ierr;
702a4c71feSBarry Smith   PetscInt       i,j,jj,jB,irow,m=A->rmap.n,*ajj,*bjj,countA,countB,colA_start,jcol;
71c1490034SHong Zhang   PetscInt       *row,*col;
72397b6df1SKris Buschelman   PetscScalar    *av, *bv,*val;
73f0c56d0fSKris Buschelman   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
74397b6df1SKris Buschelman 
75397b6df1SKris Buschelman   PetscFunctionBegin;
76397b6df1SKris Buschelman   if (mumps->isAIJ){
77397b6df1SKris Buschelman     Mat_MPIAIJ    *mat =  (Mat_MPIAIJ*)A->data;
78397b6df1SKris Buschelman     Mat_SeqAIJ    *aa=(Mat_SeqAIJ*)(mat->A)->data;
79397b6df1SKris Buschelman     Mat_SeqAIJ    *bb=(Mat_SeqAIJ*)(mat->B)->data;
80397b6df1SKris Buschelman     nz = aa->nz + bb->nz;
812a4c71feSBarry Smith     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap.rstart;
82397b6df1SKris Buschelman     garray = mat->garray;
83397b6df1SKris Buschelman     av=aa->a; bv=bb->a;
84397b6df1SKris Buschelman 
85397b6df1SKris Buschelman   } else {
86397b6df1SKris Buschelman     Mat_MPISBAIJ  *mat =  (Mat_MPISBAIJ*)A->data;
87397b6df1SKris Buschelman     Mat_SeqSBAIJ  *aa=(Mat_SeqSBAIJ*)(mat->A)->data;
88397b6df1SKris Buschelman     Mat_SeqBAIJ    *bb=(Mat_SeqBAIJ*)(mat->B)->data;
890c0e133fSBarry Smith     if (A->rmap.bs > 1) SETERRQ1(PETSC_ERR_SUP," bs=%d is not supported yet\n", A->rmap.bs);
906c6c5352SBarry Smith     nz = aa->nz + bb->nz;
912a4c71feSBarry Smith     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap.rstart;
92397b6df1SKris Buschelman     garray = mat->garray;
93397b6df1SKris Buschelman     av=aa->a; bv=bb->a;
94397b6df1SKris Buschelman   }
95397b6df1SKris Buschelman 
96397b6df1SKris Buschelman   if (!valOnly){
977c307921SBarry Smith     ierr = PetscMalloc(nz*sizeof(PetscInt) ,&row);CHKERRQ(ierr);
987c307921SBarry Smith     ierr = PetscMalloc(nz*sizeof(PetscInt),&col);CHKERRQ(ierr);
99397b6df1SKris Buschelman     ierr = PetscMalloc(nz*sizeof(PetscScalar),&val);CHKERRQ(ierr);
100397b6df1SKris Buschelman     *r = row; *c = col; *v = val;
101397b6df1SKris Buschelman   } else {
102397b6df1SKris Buschelman     row = *r; col = *c; val = *v;
103397b6df1SKris Buschelman   }
104397b6df1SKris Buschelman   *nnz = nz;
105397b6df1SKris Buschelman 
106028e57e8SHong Zhang   jj = 0; irow = rstart;
107397b6df1SKris Buschelman   for ( i=0; i<m; i++ ) {
108397b6df1SKris Buschelman     ajj = aj + ai[i];                 /* ptr to the beginning of this row */
109397b6df1SKris Buschelman     countA = ai[i+1] - ai[i];
110397b6df1SKris Buschelman     countB = bi[i+1] - bi[i];
111397b6df1SKris Buschelman     bjj = bj + bi[i];
112397b6df1SKris Buschelman 
113397b6df1SKris Buschelman     /* get jB, the starting local col index for the 2nd B-part */
114397b6df1SKris Buschelman     colA_start = rstart + ajj[0]; /* the smallest col index for A */
11575747be1SHong Zhang     j=-1;
11675747be1SHong Zhang     do {
11775747be1SHong Zhang       j++;
11875747be1SHong Zhang       if (j == countB) break;
119397b6df1SKris Buschelman       jcol = garray[bjj[j]];
12075747be1SHong Zhang     } while (jcol < colA_start);
12175747be1SHong Zhang     jB = j;
122397b6df1SKris Buschelman 
123397b6df1SKris Buschelman     /* B-part, smaller col index */
124397b6df1SKris Buschelman     colA_start = rstart + ajj[0]; /* the smallest col index for A */
125397b6df1SKris Buschelman     for (j=0; j<jB; j++){
126397b6df1SKris Buschelman       jcol = garray[bjj[j]];
127397b6df1SKris Buschelman       if (!valOnly){
128397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = jcol + shift;
12975747be1SHong Zhang 
130397b6df1SKris Buschelman       }
131397b6df1SKris Buschelman       val[jj++] = *bv++;
132397b6df1SKris Buschelman     }
133397b6df1SKris Buschelman     /* A-part */
134397b6df1SKris Buschelman     for (j=0; j<countA; j++){
135397b6df1SKris Buschelman       if (!valOnly){
136397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
137397b6df1SKris Buschelman       }
138397b6df1SKris Buschelman       val[jj++] = *av++;
139397b6df1SKris Buschelman     }
140397b6df1SKris Buschelman     /* B-part, larger col index */
141397b6df1SKris Buschelman     for (j=jB; j<countB; j++){
142397b6df1SKris Buschelman       if (!valOnly){
143397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
144397b6df1SKris Buschelman       }
145397b6df1SKris Buschelman       val[jj++] = *bv++;
146397b6df1SKris Buschelman     }
147397b6df1SKris Buschelman     irow++;
148397b6df1SKris Buschelman   }
149397b6df1SKris Buschelman 
150397b6df1SKris Buschelman   PetscFunctionReturn(0);
151397b6df1SKris Buschelman }
152397b6df1SKris Buschelman 
153c338a77dSKris Buschelman EXTERN_C_BEGIN
154c338a77dSKris Buschelman #undef __FUNCT__
155c338a77dSKris Buschelman #define __FUNCT__ "MatConvert_MUMPS_Base"
15675179d2cSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_MUMPS_Base(Mat A,MatType type,MatReuse reuse,Mat *newmat) \
157521d7252SBarry Smith {
158dfbe8321SBarry Smith   PetscErrorCode ierr;
159c338a77dSKris Buschelman   Mat            B=*newmat;
160f0c56d0fSKris Buschelman   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
161901853e0SKris Buschelman   void           (*f)(void);
162c338a77dSKris Buschelman 
163c338a77dSKris Buschelman   PetscFunctionBegin;
164ceb03754SKris Buschelman   if (reuse == MAT_INITIAL_MATRIX) {
165c338a77dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
166c338a77dSKris Buschelman   }
167f0c56d0fSKris Buschelman   B->ops->duplicate              = mumps->MatDuplicate;
168f0c56d0fSKris Buschelman   B->ops->view                   = mumps->MatView;
169f0c56d0fSKris Buschelman   B->ops->assemblyend            = mumps->MatAssemblyEnd;
170f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic       = mumps->MatLUFactorSymbolic;
171f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic = mumps->MatCholeskyFactorSymbolic;
172f0c56d0fSKris Buschelman   B->ops->destroy                = mumps->MatDestroy;
173901853e0SKris Buschelman 
174f68b968cSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr);
175901853e0SKris Buschelman   if (f) {
176f68b968cSBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C","",(PetscVoidFunction)mumps->MatPreallocate);CHKERRQ(ierr);
177901853e0SKris Buschelman   }
178901853e0SKris Buschelman 
179901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_aijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
180901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_aijmumps_seqaij_C","",PETSC_NULL);CHKERRQ(ierr);
181901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpiaij_aijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
182901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_aijmumps_mpiaij_C","",PETSC_NULL);CHKERRQ(ierr);
1832895b8caSSatish Balay   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
1842895b8caSSatish Balay   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C","",PETSC_NULL);CHKERRQ(ierr);
185901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
186901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C","",PETSC_NULL);CHKERRQ(ierr);
187901853e0SKris Buschelman 
188901853e0SKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,type);CHKERRQ(ierr);
189c338a77dSKris Buschelman   *newmat = B;
190c338a77dSKris Buschelman   PetscFunctionReturn(0);
191c338a77dSKris Buschelman }
192c338a77dSKris Buschelman EXTERN_C_END
193c338a77dSKris Buschelman 
194397b6df1SKris Buschelman #undef __FUNCT__
1953924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS"
196dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A)
197dfbe8321SBarry Smith {
198f0c56d0fSKris Buschelman   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
199dfbe8321SBarry Smith   PetscErrorCode ierr;
200c1490034SHong Zhang   PetscMPIInt    size=lu->size;
2016849ba73SBarry Smith   PetscErrorCode (*specialdestroy)(Mat);
202397b6df1SKris Buschelman   PetscFunctionBegin;
203397b6df1SKris Buschelman   if (lu->CleanUpMUMPS) {
204397b6df1SKris Buschelman     /* Terminate instance, deallocate memories */
205*329ec9b3SHong Zhang     if (size > 1){
206*329ec9b3SHong Zhang       ierr = PetscFree(lu->id.sol_loc);CHKERRQ(ierr);
207*329ec9b3SHong Zhang       ierr = VecScatterDestroy(lu->scat_rhs);CHKERRQ(ierr);
208*329ec9b3SHong Zhang       ierr = VecDestroy(lu->b_seq);CHKERRQ(ierr);
209*329ec9b3SHong Zhang       ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);
210*329ec9b3SHong Zhang       ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);
211*329ec9b3SHong Zhang       ierr = PetscFree(lu->val);CHKERRQ(ierr);
212*329ec9b3SHong Zhang     }
213397b6df1SKris Buschelman     lu->id.job=JOB_END;
214397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
215397b6df1SKris Buschelman     zmumps_c(&lu->id);
216397b6df1SKris Buschelman #else
217397b6df1SKris Buschelman     dmumps_c(&lu->id);
218397b6df1SKris Buschelman #endif
219c338a77dSKris Buschelman     ierr = PetscFree(lu->irn);CHKERRQ(ierr);
220c338a77dSKris Buschelman     ierr = PetscFree(lu->jcn);CHKERRQ(ierr);
221397b6df1SKris Buschelman     ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr);
222397b6df1SKris Buschelman   }
223a39386dcSKris Buschelman   specialdestroy = lu->specialdestroy;
224a39386dcSKris Buschelman   ierr = (*specialdestroy)(A);CHKERRQ(ierr);
225c338a77dSKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
226397b6df1SKris Buschelman   PetscFunctionReturn(0);
227397b6df1SKris Buschelman }
228397b6df1SKris Buschelman 
229397b6df1SKris Buschelman #undef __FUNCT__
230a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_AIJMUMPS"
231dfbe8321SBarry Smith PetscErrorCode MatDestroy_AIJMUMPS(Mat A)
232dfbe8321SBarry Smith {
2336849ba73SBarry Smith   PetscErrorCode ierr;
234*329ec9b3SHong Zhang   PetscMPIInt    size;
235a39386dcSKris Buschelman 
236a39386dcSKris Buschelman   PetscFunctionBegin;
237a39386dcSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
238a39386dcSKris Buschelman   if (size==1) {
239ceb03754SKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATSEQAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
240a39386dcSKris Buschelman   } else {
241ceb03754SKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATMPIAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
242a39386dcSKris Buschelman   }
243a39386dcSKris Buschelman   PetscFunctionReturn(0);
244a39386dcSKris Buschelman }
245a39386dcSKris Buschelman 
246a39386dcSKris Buschelman #undef __FUNCT__
247a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_SBAIJMUMPS"
248dfbe8321SBarry Smith PetscErrorCode MatDestroy_SBAIJMUMPS(Mat A)
249dfbe8321SBarry Smith {
2506849ba73SBarry Smith   PetscErrorCode ierr;
251*329ec9b3SHong Zhang   PetscMPIInt    size;
252a39386dcSKris Buschelman 
253a39386dcSKris Buschelman   PetscFunctionBegin;
254a39386dcSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
255a39386dcSKris Buschelman   if (size==1) {
256ceb03754SKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATSEQSBAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
257a39386dcSKris Buschelman   } else {
258ceb03754SKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATMPISBAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
259a39386dcSKris Buschelman   }
260a39386dcSKris Buschelman   PetscFunctionReturn(0);
261a39386dcSKris Buschelman }
262a39386dcSKris Buschelman 
263a39386dcSKris Buschelman #undef __FUNCT__
264c338a77dSKris Buschelman #define __FUNCT__ "MatFactorInfo_MUMPS"
265dfbe8321SBarry Smith PetscErrorCode MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) {
266f0c56d0fSKris Buschelman   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
267dfbe8321SBarry Smith   PetscErrorCode ierr;
268397b6df1SKris Buschelman 
269397b6df1SKris Buschelman   PetscFunctionBegin;
270c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
271c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                  %d \n",lu->id.sym);CHKERRQ(ierr);
272c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):           %d \n",lu->id.par);CHKERRQ(ierr);
273c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):       %d \n",lu->id.ICNTL(4));CHKERRQ(ierr);
274c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):        %d \n",lu->id.ICNTL(5));CHKERRQ(ierr);
275c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):       %d \n",lu->id.ICNTL(6));CHKERRQ(ierr);
276c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (matrix ordering):         %d \n",lu->id.ICNTL(7));CHKERRQ(ierr);
277c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(9) (A/A^T x=b is solved):     %d \n",lu->id.ICNTL(9));CHKERRQ(ierr);
278c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr);
279c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):         %d \n",lu->id.ICNTL(11));CHKERRQ(ierr);
280958c9bccSBarry Smith   if (!lu->myid && lu->id.ICNTL(11)>0) {
281c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(4) (inf norm of input mat):        %g\n",lu->id.RINFOG(4));CHKERRQ(ierr);
282c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(5) (inf norm of solution):         %g\n",lu->id.RINFOG(5));CHKERRQ(ierr);
283c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(6) (inf norm of residual):         %g\n",lu->id.RINFOG(6));CHKERRQ(ierr);
284c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",lu->id.RINFOG(7),lu->id.RINFOG(8));CHKERRQ(ierr);
285c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(9) (error estimate):               %g \n",lu->id.RINFOG(9));CHKERRQ(ierr);
286c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",lu->id.RINFOG(10),lu->id.RINFOG(11));CHKERRQ(ierr);
287c338a77dSKris Buschelman 
288c338a77dSKris Buschelman   }
289c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",lu->id.ICNTL(12));CHKERRQ(ierr);
290c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",lu->id.ICNTL(13));CHKERRQ(ierr);
291adc1d99fSHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr);
292c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(15) (efficiency control):                         %d \n",lu->id.ICNTL(15));CHKERRQ(ierr);
293c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",lu->id.ICNTL(18));CHKERRQ(ierr);
294c338a77dSKris Buschelman 
295c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",lu->id.CNTL(1));CHKERRQ(ierr);
296c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr);
297c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",lu->id.CNTL(3));CHKERRQ(ierr);
29857f0c58bSHong Zhang 
29957f0c58bSHong Zhang   /* infomation local to each processor */
300958c9bccSBarry Smith   if (!lu->myid) ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);
30157f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr);
30257f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
303958c9bccSBarry Smith   if (!lu->myid) ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);
30457f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr);
30557f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
306958c9bccSBarry Smith   if (!lu->myid) ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);
30757f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr);
30857f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
309adc1d99fSHong Zhang 
310958c9bccSBarry Smith   if (!lu->myid){ /* information from the host */
311adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr);
312adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr);
313adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr);
314adc1d99fSHong Zhang 
315adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr);
316adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr);
317adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr);
318adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr);
319adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr);
320adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr);
321adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(9) (total real space store the matrix factors after analysis): %d \n",lu->id.INFOG(9));CHKERRQ(ierr);
322adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after analysis): %d \n",lu->id.INFOG(10));CHKERRQ(ierr);
323adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix): %d \n",lu->id.INFOG(11));CHKERRQ(ierr);
324adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr);
325adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr);
326adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr);
327adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr);
328adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(16) (estimated size (in million of bytes) of all MUMPS internal data for factorization after analysis: value on the most memory consuming processor): %d \n",lu->id.INFOG(16));CHKERRQ(ierr);
329adc1d99fSHong 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);
330adc1d99fSHong 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);
331adc1d99fSHong 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);
332adc1d99fSHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr);
333adc1d99fSHong Zhang   }
334adc1d99fSHong Zhang 
335397b6df1SKris Buschelman   PetscFunctionReturn(0);
336397b6df1SKris Buschelman }
337397b6df1SKris Buschelman 
338397b6df1SKris Buschelman #undef __FUNCT__
339c1490034SHong Zhang #define __FUNCT__ "MatView_MUMPS"
340c1490034SHong Zhang PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer) {
341dfbe8321SBarry Smith   PetscErrorCode    ierr;
34232077d6dSBarry Smith   PetscTruth        iascii;
343397b6df1SKris Buschelman   PetscViewerFormat format;
344f0c56d0fSKris Buschelman   Mat_MUMPS         *mumps=(Mat_MUMPS*)(A->spptr);
345397b6df1SKris Buschelman 
346397b6df1SKris Buschelman   PetscFunctionBegin;
347397b6df1SKris Buschelman   ierr = (*mumps->MatView)(A,viewer);CHKERRQ(ierr);
348397b6df1SKris Buschelman 
34932077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
35032077d6dSBarry Smith   if (iascii) {
351397b6df1SKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
352397b6df1SKris Buschelman     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
353397b6df1SKris Buschelman       ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr);
354397b6df1SKris Buschelman     }
355397b6df1SKris Buschelman   }
356397b6df1SKris Buschelman   PetscFunctionReturn(0);
357397b6df1SKris Buschelman }
358397b6df1SKris Buschelman 
359397b6df1SKris Buschelman #undef __FUNCT__
360f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_AIJMUMPS"
361dfbe8321SBarry Smith PetscErrorCode MatSolve_AIJMUMPS(Mat A,Vec b,Vec x) {
362f0c56d0fSKris Buschelman   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
363d54de34fSKris Buschelman   PetscScalar    *array;
364397b6df1SKris Buschelman   Vec            x_seq;
365*329ec9b3SHong Zhang   IS             is_iden,is_petsc;
366*329ec9b3SHong Zhang   VecScatter     scat_rhs=lu->scat_rhs,scat_sol=lu->scat_sol;
367dfbe8321SBarry Smith   PetscErrorCode ierr;
368*329ec9b3SHong Zhang   PetscInt       i;
369397b6df1SKris Buschelman 
370397b6df1SKris Buschelman   PetscFunctionBegin;
371*329ec9b3SHong Zhang   lu->id.nrhs = 1;
372*329ec9b3SHong Zhang   x_seq = lu->b_seq;
373397b6df1SKris Buschelman   if (lu->size > 1){
374*329ec9b3SHong Zhang     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
375*329ec9b3SHong Zhang     ierr = VecScatterBegin(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat_rhs);CHKERRQ(ierr);
376*329ec9b3SHong Zhang     ierr = VecScatterEnd(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat_rhs);CHKERRQ(ierr);
377397b6df1SKris Buschelman     if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);}
378397b6df1SKris Buschelman   } else {  /* size == 1 */
379397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
380397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
381397b6df1SKris Buschelman   }
382397b6df1SKris Buschelman   if (!lu->myid) { /* define rhs on the host */
383397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
384397b6df1SKris Buschelman     lu->id.rhs = (mumps_double_complex*)array;
385397b6df1SKris Buschelman #else
386397b6df1SKris Buschelman     lu->id.rhs = array;
387397b6df1SKris Buschelman #endif
388397b6df1SKris Buschelman   }
389*329ec9b3SHong Zhang   if (lu->size == 1){
390*329ec9b3SHong Zhang     ierr = VecRestoreArray(x,&array);CHKERRQ(ierr);
391*329ec9b3SHong Zhang   } else if (!lu->myid){
392*329ec9b3SHong Zhang     ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr);
393*329ec9b3SHong Zhang   }
394*329ec9b3SHong Zhang 
395*329ec9b3SHong Zhang   if (lu->size > 1){
396*329ec9b3SHong Zhang     /* distributed solution */
397*329ec9b3SHong Zhang     lu->id.ICNTL(21) = 1;
398*329ec9b3SHong Zhang     if (!lu->nSolve){
399*329ec9b3SHong Zhang       /* Create x_seq=sol_loc for repeated use */
400*329ec9b3SHong Zhang       PetscInt    lsol_loc;
401*329ec9b3SHong Zhang       PetscScalar *sol_loc;
402*329ec9b3SHong Zhang       lsol_loc = lu->id.INFO(23); /* length of sol_loc */
403*329ec9b3SHong Zhang       ierr = PetscMalloc((1+lsol_loc)*(sizeof(PetscScalar)+sizeof(PetscInt)),&sol_loc);CHKERRQ(ierr);
404*329ec9b3SHong Zhang       lu->id.isol_loc = (PetscInt *)(sol_loc + lsol_loc);
405*329ec9b3SHong Zhang       lu->id.lsol_loc = lsol_loc;
406*329ec9b3SHong Zhang       lu->id.sol_loc  = sol_loc;
407*329ec9b3SHong Zhang       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,lsol_loc,sol_loc,&lu->x_seq);CHKERRQ(ierr);
408*329ec9b3SHong Zhang     }
409*329ec9b3SHong Zhang   }
410397b6df1SKris Buschelman 
411397b6df1SKris Buschelman   /* solve phase */
412*329ec9b3SHong Zhang   /*-------------*/
413397b6df1SKris Buschelman   lu->id.job = 3;
414397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
415397b6df1SKris Buschelman   zmumps_c(&lu->id);
416397b6df1SKris Buschelman #else
417397b6df1SKris Buschelman   dmumps_c(&lu->id);
418397b6df1SKris Buschelman #endif
419397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
42079a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1));
421397b6df1SKris Buschelman   }
422397b6df1SKris Buschelman 
423*329ec9b3SHong Zhang   if (lu->size > 1) { /* convert mumps distributed solution to petsc mpi x */
424*329ec9b3SHong Zhang     if (!lu->nSolve){ /* create scatter scat_sol */
425*329ec9b3SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,lu->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */
426*329ec9b3SHong Zhang       for (i=0; i<lu->id.lsol_loc; i++){
427*329ec9b3SHong Zhang         lu->id.isol_loc[i] -= 1; /* change Fortran style to C style */
428397b6df1SKris Buschelman       }
429*329ec9b3SHong Zhang       ierr = ISCreateGeneral(PETSC_COMM_SELF,lu->id.lsol_loc,lu->id.isol_loc,&is_petsc);CHKERRQ(ierr);  /* to */
430*329ec9b3SHong Zhang       ierr = VecScatterCreate(lu->x_seq,is_iden,x,is_petsc,&lu->scat_sol);CHKERRQ(ierr);
431*329ec9b3SHong Zhang       ierr = ISDestroy(is_iden);CHKERRQ(ierr);
432*329ec9b3SHong Zhang       ierr = ISDestroy(is_petsc);CHKERRQ(ierr);
433397b6df1SKris Buschelman     }
434*329ec9b3SHong Zhang     ierr = VecScatterBegin(lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD,lu->scat_sol);CHKERRQ(ierr);
435*329ec9b3SHong Zhang     ierr = VecScatterEnd(lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD,lu->scat_sol);CHKERRQ(ierr);
436*329ec9b3SHong Zhang   }
437*329ec9b3SHong Zhang   lu->nSolve++;
438397b6df1SKris Buschelman   PetscFunctionReturn(0);
439397b6df1SKris Buschelman }
440397b6df1SKris Buschelman 
441a58c3f20SHong Zhang /*
442a58c3f20SHong Zhang   input:
443a58c3f20SHong Zhang    F:        numeric factor
444a58c3f20SHong Zhang   output:
445a58c3f20SHong Zhang    nneg:     total number of negative pivots
446a58c3f20SHong Zhang    nzero:    0
447a58c3f20SHong Zhang    npos:     (global dimension of F) - nneg
448a58c3f20SHong Zhang */
449a58c3f20SHong Zhang 
450a58c3f20SHong Zhang #undef __FUNCT__
451a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
452dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
453a58c3f20SHong Zhang {
454a58c3f20SHong Zhang   Mat_MUMPS      *lu =(Mat_MUMPS*)F->spptr;
455dfbe8321SBarry Smith   PetscErrorCode ierr;
456c1490034SHong Zhang   PetscMPIInt    size;
457a58c3f20SHong Zhang 
458a58c3f20SHong Zhang   PetscFunctionBegin;
459bcb30aebSHong Zhang   ierr = MPI_Comm_size(F->comm,&size);CHKERRQ(ierr);
460bcb30aebSHong 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 */
461bcb30aebSHong Zhang   if (size > 1 && lu->id.ICNTL(13) != 1){
46279a5c55eSBarry Smith     SETERRQ1(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));
463bcb30aebSHong Zhang   }
464a58c3f20SHong Zhang   if (nneg){
465a58c3f20SHong Zhang     if (!lu->myid){
466a58c3f20SHong Zhang       *nneg = lu->id.INFOG(12);
467a58c3f20SHong Zhang     }
468bcb30aebSHong Zhang     ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr);
469a58c3f20SHong Zhang   }
470a58c3f20SHong Zhang   if (nzero) *nzero = 0;
4712a4c71feSBarry Smith   if (npos)  *npos  = F->rmap.N - (*nneg);
472a58c3f20SHong Zhang   PetscFunctionReturn(0);
473a58c3f20SHong Zhang }
474a58c3f20SHong Zhang 
475397b6df1SKris Buschelman #undef __FUNCT__
47619facb7aSBarry Smith #define __FUNCT__ "MatFactorNumeric_AIJMUMPS"
477af281ebdSHong Zhang PetscErrorCode MatFactorNumeric_AIJMUMPS(Mat A,MatFactorInfo *info,Mat *F)
478af281ebdSHong Zhang {
479f0c56d0fSKris Buschelman   Mat_MUMPS      *lu =(Mat_MUMPS*)(*F)->spptr;
480f0c56d0fSKris Buschelman   Mat_MUMPS      *lua=(Mat_MUMPS*)(A)->spptr;
4816849ba73SBarry Smith   PetscErrorCode ierr;
4822a4c71feSBarry Smith   PetscInt       rnz,nnz,nz=0,i,M=A->rmap.N,*ai,*aj,icntl;
483397b6df1SKris Buschelman   PetscTruth     valOnly,flg;
484e09efc27SHong Zhang   Mat            F_diag;
485397b6df1SKris Buschelman 
486397b6df1SKris Buschelman   PetscFunctionBegin;
487397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
488f0c56d0fSKris Buschelman     (*F)->ops->solve    = MatSolve_AIJMUMPS;
489397b6df1SKris Buschelman 
490397b6df1SKris Buschelman     /* Initialize a MUMPS instance */
491397b6df1SKris Buschelman     ierr = MPI_Comm_rank(A->comm, &lu->myid);
492397b6df1SKris Buschelman     ierr = MPI_Comm_size(A->comm,&lu->size);CHKERRQ(ierr);
49375747be1SHong Zhang     lua->myid = lu->myid; lua->size = lu->size;
494397b6df1SKris Buschelman     lu->id.job = JOB_INIT;
495397b6df1SKris Buschelman     ierr = MPI_Comm_dup(A->comm,&(lu->comm_mumps));CHKERRQ(ierr);
496a0e2756fSSatish Balay     ierr = MPICCommToFortranComm(lu->comm_mumps,&(lu->id.comm_fortran));CHKERRQ(ierr);
497397b6df1SKris Buschelman 
498397b6df1SKris Buschelman     /* Set mumps options */
499397b6df1SKris Buschelman     ierr = PetscOptionsBegin(A->comm,A->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
500397b6df1SKris Buschelman     lu->id.par=1;  /* host participates factorizaton and solve */
501397b6df1SKris Buschelman     lu->id.sym=lu->sym;
502397b6df1SKris Buschelman     if (lu->sym == 2){
503397b6df1SKris Buschelman       ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr);
504397b6df1SKris Buschelman       if (flg && icntl == 1) lu->id.sym=icntl;  /* matrix is spd */
505397b6df1SKris Buschelman     }
506397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
507397b6df1SKris Buschelman     zmumps_c(&lu->id);
508397b6df1SKris Buschelman #else
509397b6df1SKris Buschelman     dmumps_c(&lu->id);
510397b6df1SKris Buschelman #endif
511397b6df1SKris Buschelman 
512397b6df1SKris Buschelman     if (lu->size == 1){
513397b6df1SKris Buschelman       lu->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
514397b6df1SKris Buschelman     } else {
515397b6df1SKris Buschelman       lu->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
516397b6df1SKris Buschelman     }
517397b6df1SKris Buschelman 
518397b6df1SKris Buschelman     icntl=-1;
519397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
52019facb7aSBarry Smith     if ((flg && icntl > 0) || PetscLogPrintInfo) {
521397b6df1SKris Buschelman       lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */
522397b6df1SKris Buschelman     } else { /* no output */
523397b6df1SKris Buschelman       lu->id.ICNTL(1) = 0;  /* error message, default= 6 */
524397b6df1SKris Buschelman       lu->id.ICNTL(2) = -1; /* output stream for diagnostic printing, statistics, and warning. default=0 */
525397b6df1SKris Buschelman       lu->id.ICNTL(3) = -1; /* output stream for global information, default=6 */
526397b6df1SKris Buschelman       lu->id.ICNTL(4) = 0;  /* level of printing, 0,1,2,3,4, default=2 */
527397b6df1SKris Buschelman     }
528397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_6","ICNTL(6): matrix prescaling (0 to 7)","None",lu->id.ICNTL(6),&lu->id.ICNTL(6),PETSC_NULL);CHKERRQ(ierr);
529397b6df1SKris Buschelman     icntl=-1;
530397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
531397b6df1SKris Buschelman     if (flg) {
532397b6df1SKris Buschelman       if (icntl== 1){
533397b6df1SKris Buschelman         SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n");
534397b6df1SKris Buschelman       } else {
535397b6df1SKris Buschelman         lu->id.ICNTL(7) = icntl;
536397b6df1SKris Buschelman       }
537397b6df1SKris Buschelman     }
538397b6df1SKris Buschelman     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);
539397b6df1SKris Buschelman     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);
54094b7f48cSBarry Smith     ierr = PetscOptionsInt("-mat_mumps_icntl_11","ICNTL(11): error analysis, a positive value returns statistics (by -ksp_view)","None",lu->id.ICNTL(11),&lu->id.ICNTL(11),PETSC_NULL);CHKERRQ(ierr);
541397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr);
542397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr);
543adc1d99fSHong Zhang     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);
544397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_15","ICNTL(15): efficiency control","None",lu->id.ICNTL(15),&lu->id.ICNTL(15),PETSC_NULL);CHKERRQ(ierr);
545397b6df1SKris Buschelman 
546397b6df1SKris Buschelman     ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",lu->id.CNTL(1),&lu->id.CNTL(1),PETSC_NULL);CHKERRQ(ierr);
547397b6df1SKris Buschelman     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);
548397b6df1SKris Buschelman     ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",lu->id.CNTL(3),&lu->id.CNTL(3),PETSC_NULL);CHKERRQ(ierr);
549397b6df1SKris Buschelman     PetscOptionsEnd();
550397b6df1SKris Buschelman   }
551397b6df1SKris Buschelman 
552397b6df1SKris Buschelman   /* define matrix A */
553397b6df1SKris Buschelman   switch (lu->id.ICNTL(18)){
554397b6df1SKris Buschelman   case 0:  /* centralized assembled matrix input (size=1) */
555397b6df1SKris Buschelman     if (!lu->myid) {
556c36ead0aSKris Buschelman       if (lua->isAIJ){
557397b6df1SKris Buschelman         Mat_SeqAIJ   *aa = (Mat_SeqAIJ*)A->data;
558397b6df1SKris Buschelman         nz               = aa->nz;
559397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
560397b6df1SKris Buschelman       } else {
561397b6df1SKris Buschelman         Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data;
5626c6c5352SBarry Smith         nz                  =  aa->nz;
563397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
564397b6df1SKris Buschelman       }
565397b6df1SKris Buschelman       if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */
5667c307921SBarry Smith         ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->irn);CHKERRQ(ierr);
5677c307921SBarry Smith         ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->jcn);CHKERRQ(ierr);
568397b6df1SKris Buschelman         nz = 0;
569397b6df1SKris Buschelman         for (i=0; i<M; i++){
570397b6df1SKris Buschelman           rnz = ai[i+1] - ai[i];
571397b6df1SKris Buschelman           while (rnz--) {  /* Fortran row/col index! */
572397b6df1SKris Buschelman             lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++;
573397b6df1SKris Buschelman           }
574397b6df1SKris Buschelman         }
575397b6df1SKris Buschelman       }
576397b6df1SKris Buschelman     }
577397b6df1SKris Buschelman     break;
578397b6df1SKris Buschelman   case 3:  /* distributed assembled matrix input (size>1) */
579397b6df1SKris Buschelman     if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
580397b6df1SKris Buschelman       valOnly = PETSC_FALSE;
581397b6df1SKris Buschelman     } else {
582397b6df1SKris Buschelman       valOnly = PETSC_TRUE; /* only update mat values, not row and col index */
583397b6df1SKris Buschelman     }
584397b6df1SKris Buschelman     ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
585397b6df1SKris Buschelman     break;
586397b6df1SKris Buschelman   default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS.");
587397b6df1SKris Buschelman   }
588397b6df1SKris Buschelman 
589397b6df1SKris Buschelman   /* analysis phase */
590*329ec9b3SHong Zhang   /*----------------*/
591397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
592*329ec9b3SHong Zhang     lu->id.job = 1;
593*329ec9b3SHong Zhang 
594397b6df1SKris Buschelman     lu->id.n = M;
595397b6df1SKris Buschelman     switch (lu->id.ICNTL(18)){
596397b6df1SKris Buschelman     case 0:  /* centralized assembled matrix input */
597397b6df1SKris Buschelman       if (!lu->myid) {
598397b6df1SKris Buschelman         lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
599397b6df1SKris Buschelman         if (lu->id.ICNTL(6)>1){
600397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
601397b6df1SKris Buschelman           lu->id.a = (mumps_double_complex*)lu->val;
602397b6df1SKris Buschelman #else
603397b6df1SKris Buschelman           lu->id.a = lu->val;
604397b6df1SKris Buschelman #endif
605397b6df1SKris Buschelman         }
606397b6df1SKris Buschelman       }
607397b6df1SKris Buschelman       break;
608397b6df1SKris Buschelman     case 3:  /* distributed assembled matrix input (size>1) */
609397b6df1SKris Buschelman       lu->id.nz_loc = nnz;
610397b6df1SKris Buschelman       lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
611397b6df1SKris Buschelman       if (lu->id.ICNTL(6)>1) {
612397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
613397b6df1SKris Buschelman         lu->id.a_loc = (mumps_double_complex*)lu->val;
614397b6df1SKris Buschelman #else
615397b6df1SKris Buschelman         lu->id.a_loc = lu->val;
616397b6df1SKris Buschelman #endif
617397b6df1SKris Buschelman       }
618*329ec9b3SHong Zhang       /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
619*329ec9b3SHong Zhang       IS  is_iden;
620*329ec9b3SHong Zhang       Vec b;
621*329ec9b3SHong Zhang       if (!lu->myid){
622*329ec9b3SHong Zhang         ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap.N,&lu->b_seq);CHKERRQ(ierr);
623*329ec9b3SHong Zhang         ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap.N,0,1,&is_iden);CHKERRQ(ierr);
624*329ec9b3SHong Zhang       } else {
625*329ec9b3SHong Zhang         ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr);
626*329ec9b3SHong Zhang         ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
627*329ec9b3SHong Zhang       }
628*329ec9b3SHong Zhang       ierr = VecCreate(A->comm,&b);CHKERRQ(ierr);
629*329ec9b3SHong Zhang       ierr = VecSetSizes(b,A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr);
630*329ec9b3SHong Zhang       ierr = VecSetFromOptions(b);CHKERRQ(ierr);
631*329ec9b3SHong Zhang 
632*329ec9b3SHong Zhang       ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr);
633*329ec9b3SHong Zhang       ierr = ISDestroy(is_iden);CHKERRQ(ierr);
634*329ec9b3SHong Zhang       ierr = VecDestroy(b);CHKERRQ(ierr);
635397b6df1SKris Buschelman       break;
636397b6df1SKris Buschelman     }
637397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
638397b6df1SKris Buschelman     zmumps_c(&lu->id);
639397b6df1SKris Buschelman #else
640397b6df1SKris Buschelman     dmumps_c(&lu->id);
641397b6df1SKris Buschelman #endif
642397b6df1SKris Buschelman     if (lu->id.INFOG(1) < 0) {
64379a5c55eSBarry Smith       SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
644397b6df1SKris Buschelman     }
645397b6df1SKris Buschelman   }
646397b6df1SKris Buschelman 
647397b6df1SKris Buschelman   /* numerical factorization phase */
648*329ec9b3SHong Zhang   /*-------------------------------*/
649*329ec9b3SHong Zhang   lu->id.job = 2;
650958c9bccSBarry Smith   if(!lu->id.ICNTL(18)) {
651a7aca84bSHong Zhang     if (!lu->myid) {
652397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
653397b6df1SKris Buschelman       lu->id.a = (mumps_double_complex*)lu->val;
654397b6df1SKris Buschelman #else
655397b6df1SKris Buschelman       lu->id.a = lu->val;
656397b6df1SKris Buschelman #endif
657397b6df1SKris Buschelman     }
658397b6df1SKris Buschelman   } else {
659397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
660397b6df1SKris Buschelman     lu->id.a_loc = (mumps_double_complex*)lu->val;
661397b6df1SKris Buschelman #else
662397b6df1SKris Buschelman     lu->id.a_loc = lu->val;
663397b6df1SKris Buschelman #endif
664397b6df1SKris Buschelman   }
665397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
666397b6df1SKris Buschelman   zmumps_c(&lu->id);
667397b6df1SKris Buschelman #else
668397b6df1SKris Buschelman   dmumps_c(&lu->id);
669397b6df1SKris Buschelman #endif
670397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
67119facb7aSBarry Smith     if (lu->id.INFO(1) == -13) {
67219facb7aSBarry Smith       SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",lu->id.INFO(2));
67319facb7aSBarry Smith     } else {
67479a5c55eSBarry Smith       SETERRQ2(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));
675397b6df1SKris Buschelman     }
67619facb7aSBarry Smith   }
677397b6df1SKris Buschelman 
67819facb7aSBarry Smith   if (!lu->myid && lu->id.ICNTL(16) > 0){
67979a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_LIB,"  lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16));
680397b6df1SKris Buschelman   }
681397b6df1SKris Buschelman 
6828ada1bb4SHong Zhang   if (lu->size > 1){
683e09efc27SHong Zhang     if ((*F)->factor == FACTOR_LU){
684e09efc27SHong Zhang       F_diag = ((Mat_MPIAIJ *)(*F)->data)->A;
685e09efc27SHong Zhang     } else {
686e09efc27SHong Zhang       F_diag = ((Mat_MPISBAIJ *)(*F)->data)->A;
687e09efc27SHong Zhang     }
688e09efc27SHong Zhang     F_diag->assembled = PETSC_TRUE;
689*329ec9b3SHong Zhang     if (lu->nSolve){
690*329ec9b3SHong Zhang       ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);
691*329ec9b3SHong Zhang       ierr = PetscFree(lu->id.sol_loc);CHKERRQ(ierr);
692*329ec9b3SHong Zhang       ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);
693*329ec9b3SHong Zhang     }
6948ada1bb4SHong Zhang   }
695397b6df1SKris Buschelman   (*F)->assembled   = PETSC_TRUE;
696397b6df1SKris Buschelman   lu->matstruc      = SAME_NONZERO_PATTERN;
697ace87b0dSHong Zhang   lu->CleanUpMUMPS  = PETSC_TRUE;
698*329ec9b3SHong Zhang   lu->nSolve        = 0;
699397b6df1SKris Buschelman   PetscFunctionReturn(0);
700397b6df1SKris Buschelman }
701397b6df1SKris Buschelman 
702397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */
703397b6df1SKris Buschelman #undef __FUNCT__
704f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
705dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) {
706397b6df1SKris Buschelman   Mat            B;
707f0c56d0fSKris Buschelman   Mat_MUMPS      *lu;
708dfbe8321SBarry Smith   PetscErrorCode ierr;
709397b6df1SKris Buschelman 
710397b6df1SKris Buschelman   PetscFunctionBegin;
711397b6df1SKris Buschelman   /* Create the factorization matrix */
712f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&B);CHKERRQ(ierr);
7132a4c71feSBarry Smith   ierr = MatSetSizes(B,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);CHKERRQ(ierr);
714be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
715397b6df1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
716397b6df1SKris Buschelman   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
717397b6df1SKris Buschelman 
718f0c56d0fSKris Buschelman   B->ops->lufactornumeric = MatFactorNumeric_AIJMUMPS;
719397b6df1SKris Buschelman   B->factor               = FACTOR_LU;
720f0c56d0fSKris Buschelman   lu                      = (Mat_MUMPS*)B->spptr;
721397b6df1SKris Buschelman   lu->sym                 = 0;
722397b6df1SKris Buschelman   lu->matstruc            = DIFFERENT_NONZERO_PATTERN;
723397b6df1SKris Buschelman 
724397b6df1SKris Buschelman   *F = B;
725397b6df1SKris Buschelman   PetscFunctionReturn(0);
726397b6df1SKris Buschelman }
727397b6df1SKris Buschelman 
728397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */
729397b6df1SKris Buschelman #undef __FUNCT__
730f0c56d0fSKris Buschelman #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS"
731dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat A,IS r,MatFactorInfo *info,Mat *F) {
732397b6df1SKris Buschelman   Mat            B;
733f0c56d0fSKris Buschelman   Mat_MUMPS      *lu;
734dfbe8321SBarry Smith   PetscErrorCode ierr;
735397b6df1SKris Buschelman 
736397b6df1SKris Buschelman   PetscFunctionBegin;
737397b6df1SKris Buschelman   /* Create the factorization matrix */
738f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&B);CHKERRQ(ierr);
7392a4c71feSBarry Smith   ierr = MatSetSizes(B,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);CHKERRQ(ierr);
740be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
741efc670deSHong Zhang   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
742efc670deSHong Zhang   ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
743397b6df1SKris Buschelman 
744f0c56d0fSKris Buschelman   B->ops->choleskyfactornumeric = MatFactorNumeric_AIJMUMPS;
745a58c3f20SHong Zhang   B->ops->getinertia            = MatGetInertia_SBAIJMUMPS;
746397b6df1SKris Buschelman   B->factor                     = FACTOR_CHOLESKY;
747f0c56d0fSKris Buschelman   lu                            = (Mat_MUMPS*)B->spptr;
748397b6df1SKris Buschelman   lu->sym                       = 2;
749397b6df1SKris Buschelman   lu->matstruc                  = DIFFERENT_NONZERO_PATTERN;
750397b6df1SKris Buschelman 
751397b6df1SKris Buschelman   *F = B;
752397b6df1SKris Buschelman   PetscFunctionReturn(0);
753397b6df1SKris Buschelman }
754397b6df1SKris Buschelman 
755397b6df1SKris Buschelman #undef __FUNCT__
756f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_AIJMUMPS"
757dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_AIJMUMPS(Mat A,MatAssemblyType mode) {
758dfbe8321SBarry Smith   PetscErrorCode ierr;
759f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
760c338a77dSKris Buschelman 
761397b6df1SKris Buschelman   PetscFunctionBegin;
762c338a77dSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
763f0c56d0fSKris Buschelman 
764c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
765c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
766f0c56d0fSKris Buschelman   A->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
767397b6df1SKris Buschelman   PetscFunctionReturn(0);
768397b6df1SKris Buschelman }
769397b6df1SKris Buschelman 
770c338a77dSKris Buschelman EXTERN_C_BEGIN
771c338a77dSKris Buschelman #undef __FUNCT__
772f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_AIJ_AIJMUMPS"
77375179d2cSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_AIJ_AIJMUMPS(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
774dfbe8321SBarry Smith {
775dfbe8321SBarry Smith   PetscErrorCode ierr;
776521d7252SBarry Smith   PetscMPIInt    size;
777c338a77dSKris Buschelman   MPI_Comm       comm;
778c338a77dSKris Buschelman   Mat            B=*newmat;
779f0c56d0fSKris Buschelman   Mat_MUMPS      *mumps;
780397b6df1SKris Buschelman 
781397b6df1SKris Buschelman   PetscFunctionBegin;
782ceb03754SKris Buschelman   if (reuse == MAT_INITIAL_MATRIX) {
783c338a77dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
784397b6df1SKris Buschelman   }
785397b6df1SKris Buschelman 
786c338a77dSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
787f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
788c338a77dSKris Buschelman 
789f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
790c338a77dSKris Buschelman   mumps->MatView                   = A->ops->view;
791c338a77dSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
792c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
793c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
794c338a77dSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
795a39386dcSKris Buschelman   mumps->specialdestroy            = MatDestroy_AIJMUMPS;
796c338a77dSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
797f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_TRUE;
798c338a77dSKris Buschelman 
7994b68dd72SKris Buschelman   B->spptr                         = (void*)mumps;
800422e82a1SHong Zhang   B->ops->duplicate                = MatDuplicate_MUMPS;
801c1490034SHong Zhang   B->ops->view                     = MatView_MUMPS;
802f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_AIJMUMPS;
803f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
8043924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
805c338a77dSKris Buschelman 
806c338a77dSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
807c338a77dSKris Buschelman   if (size == 1) {
808c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_aijmumps_C",
809f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
810c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_seqaij_C",
811c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
812c338a77dSKris Buschelman   } else {
813c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpiaij_aijmumps_C",
814f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
815c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_mpiaij_C",
816c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
817c338a77dSKris Buschelman   }
818c338a77dSKris Buschelman 
819ae15b995SBarry Smith   ierr = PetscInfo(0,"Using MUMPS for LU factorization and solves.\n");CHKERRQ(ierr);
820c338a77dSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
821c338a77dSKris Buschelman   *newmat = B;
822397b6df1SKris Buschelman   PetscFunctionReturn(0);
823397b6df1SKris Buschelman }
824c338a77dSKris Buschelman EXTERN_C_END
825397b6df1SKris Buschelman 
82624b6179bSKris Buschelman /*MC
827fafad747SKris Buschelman   MATAIJMUMPS - MATAIJMUMPS = "aijmumps" - A matrix type providing direct solvers (LU) for distributed
82824b6179bSKris Buschelman   and sequential matrices via the external package MUMPS.
82924b6179bSKris Buschelman 
83024b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
83124b6179bSKris Buschelman   on how to declare the existence of external packages),
83224b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
83324b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATAIJMUMPS).
83424b6179bSKris Buschelman 
83524b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQAIJ.
83624b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPIAIJ.  Hence for single process communicators,
83724b6179bSKris Buschelman   MatSeqAIJSetPreallocation is supported, and similarly MatMPIAIJSetPreallocation is supported
83824b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
83928b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
84028b08bd3SKris Buschelman   conversion to or from the MATSEQAIJ or MATMPIAIJ type (depending on the communicator size)
84128b08bd3SKris Buschelman   without data copy.
84224b6179bSKris Buschelman 
84324b6179bSKris Buschelman   Options Database Keys:
8440bad9183SKris Buschelman + -mat_type aijmumps - sets the matrix type to "aijmumps" during a call to MatSetFromOptions()
84524b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
84624b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,1,2,3,4> - print level
84724b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
84824b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
84924b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
85024b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
85194b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
85224b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
85324b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
85424b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
85524b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
85624b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
85724b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
85824b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
85924b6179bSKris Buschelman 
86024b6179bSKris Buschelman   Level: beginner
86124b6179bSKris Buschelman 
86224b6179bSKris Buschelman .seealso: MATSBAIJMUMPS
86324b6179bSKris Buschelman M*/
86424b6179bSKris Buschelman 
865397b6df1SKris Buschelman EXTERN_C_BEGIN
866397b6df1SKris Buschelman #undef __FUNCT__
867f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_AIJMUMPS"
868be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_AIJMUMPS(Mat A)
869dfbe8321SBarry Smith {
870dfbe8321SBarry Smith   PetscErrorCode ierr;
871c1490034SHong Zhang   PetscMPIInt    size;
872397b6df1SKris Buschelman   MPI_Comm       comm;
873397b6df1SKris Buschelman 
874397b6df1SKris Buschelman   PetscFunctionBegin;
8755441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqAIJ or MPIAIJ */
8765441df8eSKris Buschelman   /*   and AIJMUMPS types */
8775441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATAIJMUMPS);CHKERRQ(ierr);
878397b6df1SKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
879397b6df1SKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
880397b6df1SKris Buschelman   if (size == 1) {
881397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
882397b6df1SKris Buschelman   } else {
883397b6df1SKris Buschelman     ierr   = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
88400ff2a26SHong Zhang     /*
88500ff2a26SHong Zhang     Mat A_diag = ((Mat_MPIAIJ *)A->data)->A;
886ceb03754SKris Buschelman     ierr   = MatConvert_AIJ_AIJMUMPS(A_diag,MATAIJMUMPS,MAT_REUSE_MATRIX,&A_diag);CHKERRQ(ierr);
88700ff2a26SHong Zhang     */
888397b6df1SKris Buschelman   }
889ceb03754SKris Buschelman   ierr = MatConvert_AIJ_AIJMUMPS(A,MATAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
890397b6df1SKris Buschelman   PetscFunctionReturn(0);
891397b6df1SKris Buschelman }
892397b6df1SKris Buschelman EXTERN_C_END
893397b6df1SKris Buschelman 
894f579278aSKris Buschelman #undef __FUNCT__
895f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SBAIJMUMPS"
896dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SBAIJMUMPS(Mat A,MatAssemblyType mode)
897dfbe8321SBarry Smith {
898dfbe8321SBarry Smith   PetscErrorCode ierr;
899f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
900f579278aSKris Buschelman 
901f579278aSKris Buschelman   PetscFunctionBegin;
902f579278aSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
903f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
904f0c56d0fSKris Buschelman   A->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
905f579278aSKris Buschelman   PetscFunctionReturn(0);
906f579278aSKris Buschelman }
907f579278aSKris Buschelman 
908f579278aSKris Buschelman EXTERN_C_BEGIN
909f579278aSKris Buschelman #undef __FUNCT__
9109c097c71SKris Buschelman #define __FUNCT__ "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS"
911c1490034SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation_MPISBAIJMUMPS(Mat B,PetscInt bs,PetscInt d_nz,PetscInt *d_nnz,PetscInt o_nz,PetscInt *o_nnz)
9129c097c71SKris Buschelman {
9139c097c71SKris Buschelman   Mat       A;
91402217bfdSHong Zhang   Mat_MUMPS *mumps=(Mat_MUMPS*)B->spptr;
915dfbe8321SBarry Smith   PetscErrorCode ierr;
9169c097c71SKris Buschelman 
9179c097c71SKris Buschelman   PetscFunctionBegin;
9189c097c71SKris Buschelman   /*
9199c097c71SKris Buschelman     After performing the MPISBAIJ Preallocation, we need to convert the local diagonal block matrix
9209c097c71SKris Buschelman     into MUMPS type so that the block jacobi preconditioner (for example) can use MUMPS.  I would
9219c097c71SKris Buschelman     like this to be done in the MatCreate routine, but the creation of this inner matrix requires
9229c097c71SKris Buschelman     block size info so that PETSc can determine the local size properly.  The block size info is set
9239c097c71SKris Buschelman     in the preallocation routine.
9249c097c71SKris Buschelman   */
9259c097c71SKris Buschelman   ierr = (*mumps->MatPreallocate)(B,bs,d_nz,d_nnz,o_nz,o_nnz);
9269c097c71SKris Buschelman   A    = ((Mat_MPISBAIJ *)B->data)->A;
927ceb03754SKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
9289c097c71SKris Buschelman   PetscFunctionReturn(0);
9299c097c71SKris Buschelman }
9309c097c71SKris Buschelman EXTERN_C_END
9319c097c71SKris Buschelman 
9329c097c71SKris Buschelman EXTERN_C_BEGIN
9339c097c71SKris Buschelman #undef __FUNCT__
934f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_SBAIJ_SBAIJMUMPS"
93575179d2cSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SBAIJ_SBAIJMUMPS(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
936dfbe8321SBarry Smith {
937dfbe8321SBarry Smith   PetscErrorCode ierr;
938521d7252SBarry Smith   PetscMPIInt    size;
939f579278aSKris Buschelman   MPI_Comm       comm;
940f579278aSKris Buschelman   Mat            B=*newmat;
941422e82a1SHong Zhang   Mat_MUMPS      *mumps;
9429c097c71SKris Buschelman   void           (*f)(void);
943f579278aSKris Buschelman 
944f579278aSKris Buschelman   PetscFunctionBegin;
945ceb03754SKris Buschelman   if (reuse == MAT_INITIAL_MATRIX) {
946f579278aSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
947f579278aSKris Buschelman   }
948f579278aSKris Buschelman 
949f579278aSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
950f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
951f579278aSKris Buschelman 
952f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
953f579278aSKris Buschelman   mumps->MatView                   = A->ops->view;
954f579278aSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
955f579278aSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
956f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
957f579278aSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
958a39386dcSKris Buschelman   mumps->specialdestroy            = MatDestroy_SBAIJMUMPS;
959f579278aSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
960f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_FALSE;
961f579278aSKris Buschelman 
962f579278aSKris Buschelman   B->spptr                         = (void*)mumps;
963422e82a1SHong Zhang   B->ops->duplicate                = MatDuplicate_MUMPS;
964c1490034SHong Zhang   B->ops->view                     = MatView_MUMPS;
965f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_SBAIJMUMPS;
966f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
9673924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
968f579278aSKris Buschelman 
969f579278aSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
970f579278aSKris Buschelman   if (size == 1) {
971f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C",
972f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
973f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C",
974f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
975f579278aSKris Buschelman   } else {
9769c097c71SKris Buschelman   /* I really don't like needing to know the tag: MatMPISBAIJSetPreallocation_C */
977f68b968cSBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr);
978901853e0SKris Buschelman     if (f) { /* This case should always be true when this routine is called */
9796849ba73SBarry Smith       mumps->MatPreallocate = (PetscErrorCode (*)(Mat,int,int,int*,int,int*))f;
9809c097c71SKris Buschelman       ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPISBAIJSetPreallocation_C",
9819c097c71SKris Buschelman                                                "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS",
9829c097c71SKris Buschelman                                                MatMPISBAIJSetPreallocation_MPISBAIJMUMPS);CHKERRQ(ierr);
9839c097c71SKris Buschelman     }
984f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C",
985f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
986f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C",
987f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
988f579278aSKris Buschelman   }
989f579278aSKris Buschelman 
990ae15b995SBarry Smith   ierr = PetscInfo(0,"Using MUMPS for Cholesky factorization and solves.\n");CHKERRQ(ierr);
991f579278aSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
992f579278aSKris Buschelman   *newmat = B;
993f579278aSKris Buschelman   PetscFunctionReturn(0);
994f579278aSKris Buschelman }
995f579278aSKris Buschelman EXTERN_C_END
996f579278aSKris Buschelman 
997f0c56d0fSKris Buschelman #undef __FUNCT__
998422e82a1SHong Zhang #define __FUNCT__ "MatDuplicate_MUMPS"
999dfbe8321SBarry Smith PetscErrorCode MatDuplicate_MUMPS(Mat A, MatDuplicateOption op, Mat *M) {
1000dfbe8321SBarry Smith   PetscErrorCode ierr;
10018e393735SKris Buschelman   Mat_MUMPS   *lu=(Mat_MUMPS *)A->spptr;
10028f340917SKris Buschelman 
1003f0c56d0fSKris Buschelman   PetscFunctionBegin;
10048f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
10058e393735SKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_MUMPS));CHKERRQ(ierr);
1006f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
1007f0c56d0fSKris Buschelman }
1008f0c56d0fSKris Buschelman 
100924b6179bSKris Buschelman /*MC
1010fafad747SKris Buschelman   MATSBAIJMUMPS - MATSBAIJMUMPS = "sbaijmumps" - A symmetric matrix type providing direct solvers (Cholesky) for
101124b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
101224b6179bSKris Buschelman 
101324b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
101424b6179bSKris Buschelman   on how to declare the existence of external packages),
101524b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
101624b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATSBAIJMUMPS).
101724b6179bSKris Buschelman 
101824b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQSBAIJ.
101924b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPISBAIJ.  Hence for single process communicators,
102024b6179bSKris Buschelman   MatSeqSBAIJSetPreallocation is supported, and similarly MatMPISBAIJSetPreallocation is supported
102124b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
102228b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
102328b08bd3SKris Buschelman   conversion to or from the MATSEQSBAIJ or MATMPISBAIJ type (depending on the communicator size)
102428b08bd3SKris Buschelman   without data copy.
102524b6179bSKris Buschelman 
102624b6179bSKris Buschelman   Options Database Keys:
10270bad9183SKris Buschelman + -mat_type sbaijmumps - sets the matrix type to "sbaijmumps" during a call to MatSetFromOptions()
102824b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
102924b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,...,4> - print level
103024b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
103124b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
103224b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
103324b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
103494b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
103524b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
103624b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
103724b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
103824b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
103924b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
104024b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
104124b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
104224b6179bSKris Buschelman 
104324b6179bSKris Buschelman   Level: beginner
104424b6179bSKris Buschelman 
104524b6179bSKris Buschelman .seealso: MATAIJMUMPS
104624b6179bSKris Buschelman M*/
104724b6179bSKris Buschelman 
1048397b6df1SKris Buschelman EXTERN_C_BEGIN
1049397b6df1SKris Buschelman #undef __FUNCT__
1050f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SBAIJMUMPS"
1051be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SBAIJMUMPS(Mat A)
1052dfbe8321SBarry Smith {
10536849ba73SBarry Smith   PetscErrorCode ierr;
10548ada1bb4SHong Zhang   PetscMPIInt    size;
1055397b6df1SKris Buschelman 
1056397b6df1SKris Buschelman   PetscFunctionBegin;
10575441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqSBAIJ or MPISBAIJ */
10585441df8eSKris Buschelman   /*   and SBAIJMUMPS types */
10595441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSBAIJMUMPS);CHKERRQ(ierr);
10605441df8eSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
1061397b6df1SKris Buschelman   if (size == 1) {
1062397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQSBAIJ);CHKERRQ(ierr);
1063397b6df1SKris Buschelman   } else {
1064397b6df1SKris Buschelman     ierr   = MatSetType(A,MATMPISBAIJ);CHKERRQ(ierr);
1065397b6df1SKris Buschelman   }
1066ceb03754SKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
1067397b6df1SKris Buschelman   PetscFunctionReturn(0);
1068397b6df1SKris Buschelman }
1069397b6df1SKris Buschelman EXTERN_C_END
1070