xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision bd920abb044c96993aaee00a1f0a41db581d8359)
1 #define PETSCMAT_DLL
2 
3 /*
4     Provides an interface to the MUMPS sparse solver
5 */
6 #include "src/mat/impls/aij/seq/aij.h"
7 #include "src/mat/impls/aij/mpi/mpiaij.h"
8 #include "src/mat/impls/sbaij/seq/sbaij.h"
9 #include "src/mat/impls/sbaij/mpi/mpisbaij.h"
10 
11 EXTERN_C_BEGIN
12 #if defined(PETSC_USE_COMPLEX)
13 #include "zmumps_c.h"
14 #else
15 #include "dmumps_c.h"
16 #endif
17 EXTERN_C_END
18 #define JOB_INIT -1
19 #define JOB_END -2
20 /* macros s.t. indices match MUMPS documentation */
21 #define ICNTL(I) icntl[(I)-1]
22 #define CNTL(I) cntl[(I)-1]
23 #define INFOG(I) infog[(I)-1]
24 #define INFO(I) info[(I)-1]
25 #define RINFOG(I) rinfog[(I)-1]
26 #define RINFO(I) rinfo[(I)-1]
27 
28 typedef struct {
29 #if defined(PETSC_USE_COMPLEX)
30   ZMUMPS_STRUC_C id;
31 #else
32   DMUMPS_STRUC_C id;
33 #endif
34   MatStructure   matstruc;
35   PetscMPIInt    myid,size;
36   PetscInt       *irn,*jcn,sym,nSolve;
37   PetscScalar    *val;
38   MPI_Comm       comm_mumps;
39   VecScatter     scat_rhs, scat_sol;
40   PetscTruth     isAIJ,CleanUpMUMPS;
41   Vec            b_seq,x_seq;
42   PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
43   PetscErrorCode (*MatView)(Mat,PetscViewer);
44   PetscErrorCode (*MatAssemblyEnd)(Mat,MatAssemblyType);
45   PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
46   PetscErrorCode (*MatCholeskyFactorSymbolic)(Mat,IS,MatFactorInfo*,Mat*);
47   PetscErrorCode (*MatDestroy)(Mat);
48   PetscErrorCode (*specialdestroy)(Mat);
49   PetscErrorCode (*MatPreallocate)(Mat,int,int,int*,int,int*);
50 } Mat_MUMPS;
51 
52 EXTERN PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);
53 EXTERN_C_BEGIN
54 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SBAIJ_SBAIJMUMPS(Mat,MatType,MatReuse,Mat*);
55 EXTERN_C_END
56 /* convert Petsc mpiaij matrix to triples: row[nz], col[nz], val[nz] */
57 /*
58   input:
59     A       - matrix in mpiaij or mpisbaij (bs=1) format
60     shift   - 0: C style output triple; 1: Fortran style output triple.
61     valOnly - FALSE: spaces are allocated and values are set for the triple
62               TRUE:  only the values in v array are updated
63   output:
64     nnz     - dim of r, c, and v (number of local nonzero entries of A)
65     r, c, v - row and col index, matrix values (matrix triples)
66  */
67 PetscErrorCode MatConvertToTriples(Mat A,int shift,PetscTruth valOnly,int *nnz,int **r, int **c, PetscScalar **v) {
68   PetscInt       *ai, *aj, *bi, *bj, rstart,nz, *garray;
69   PetscErrorCode ierr;
70   PetscInt       i,j,jj,jB,irow,m=A->rmap.n,*ajj,*bjj,countA,countB,colA_start,jcol;
71   PetscInt       *row,*col;
72   PetscScalar    *av, *bv,*val;
73   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
74 
75   PetscFunctionBegin;
76   if (mumps->isAIJ){
77     Mat_MPIAIJ    *mat =  (Mat_MPIAIJ*)A->data;
78     Mat_SeqAIJ    *aa=(Mat_SeqAIJ*)(mat->A)->data;
79     Mat_SeqAIJ    *bb=(Mat_SeqAIJ*)(mat->B)->data;
80     nz = aa->nz + bb->nz;
81     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap.rstart;
82     garray = mat->garray;
83     av=aa->a; bv=bb->a;
84 
85   } else {
86     Mat_MPISBAIJ  *mat =  (Mat_MPISBAIJ*)A->data;
87     Mat_SeqSBAIJ  *aa=(Mat_SeqSBAIJ*)(mat->A)->data;
88     Mat_SeqBAIJ    *bb=(Mat_SeqBAIJ*)(mat->B)->data;
89     if (A->rmap.bs > 1) SETERRQ1(PETSC_ERR_SUP," bs=%d is not supported yet\n", A->rmap.bs);
90     nz = aa->nz + bb->nz;
91     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap.rstart;
92     garray = mat->garray;
93     av=aa->a; bv=bb->a;
94   }
95 
96   if (!valOnly){
97     ierr = PetscMalloc(nz*sizeof(PetscInt) ,&row);CHKERRQ(ierr);
98     ierr = PetscMalloc(nz*sizeof(PetscInt),&col);CHKERRQ(ierr);
99     ierr = PetscMalloc(nz*sizeof(PetscScalar),&val);CHKERRQ(ierr);
100     *r = row; *c = col; *v = val;
101   } else {
102     row = *r; col = *c; val = *v;
103   }
104   *nnz = nz;
105 
106   jj = 0; irow = rstart;
107   for ( i=0; i<m; i++ ) {
108     ajj = aj + ai[i];                 /* ptr to the beginning of this row */
109     countA = ai[i+1] - ai[i];
110     countB = bi[i+1] - bi[i];
111     bjj = bj + bi[i];
112 
113     /* get jB, the starting local col index for the 2nd B-part */
114     colA_start = rstart + ajj[0]; /* the smallest col index for A */
115     j=-1;
116     do {
117       j++;
118       if (j == countB) break;
119       jcol = garray[bjj[j]];
120     } while (jcol < colA_start);
121     jB = j;
122 
123     /* B-part, smaller col index */
124     colA_start = rstart + ajj[0]; /* the smallest col index for A */
125     for (j=0; j<jB; j++){
126       jcol = garray[bjj[j]];
127       if (!valOnly){
128         row[jj] = irow + shift; col[jj] = jcol + shift;
129 
130       }
131       val[jj++] = *bv++;
132     }
133     /* A-part */
134     for (j=0; j<countA; j++){
135       if (!valOnly){
136         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
137       }
138       val[jj++] = *av++;
139     }
140     /* B-part, larger col index */
141     for (j=jB; j<countB; j++){
142       if (!valOnly){
143         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
144       }
145       val[jj++] = *bv++;
146     }
147     irow++;
148   }
149 
150   PetscFunctionReturn(0);
151 }
152 
153 EXTERN_C_BEGIN
154 #undef __FUNCT__
155 #define __FUNCT__ "MatConvert_MUMPS_Base"
156 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_MUMPS_Base(Mat A,MatType type,MatReuse reuse,Mat *newmat)
157 {
158   PetscErrorCode ierr;
159   Mat            B=*newmat;
160   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
161   void           (*f)(void);
162 
163   PetscFunctionBegin;
164   if (reuse == MAT_INITIAL_MATRIX) {
165     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
166   }
167   B->ops->duplicate              = mumps->MatDuplicate;
168   B->ops->view                   = mumps->MatView;
169   B->ops->assemblyend            = mumps->MatAssemblyEnd;
170   B->ops->lufactorsymbolic       = mumps->MatLUFactorSymbolic;
171   B->ops->choleskyfactorsymbolic = mumps->MatCholeskyFactorSymbolic;
172   B->ops->destroy                = mumps->MatDestroy;
173 
174   /* put back original composed preallocation function */
175   ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr);
176   if (f) {
177     ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C","",(PetscVoidFunction)mumps->MatPreallocate);CHKERRQ(ierr);
178   }
179   ierr     = PetscFree(mumps);CHKERRQ(ierr);
180   A->spptr = PETSC_NULL;
181 
182   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_aijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
183   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_aijmumps_seqaij_C","",PETSC_NULL);CHKERRQ(ierr);
184   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpiaij_aijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
185   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_aijmumps_mpiaij_C","",PETSC_NULL);CHKERRQ(ierr);
186   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
187   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C","",PETSC_NULL);CHKERRQ(ierr);
188   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
189   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C","",PETSC_NULL);CHKERRQ(ierr);
190 
191   ierr = PetscObjectChangeTypeName((PetscObject)B,type);CHKERRQ(ierr);
192   *newmat = B;
193   PetscFunctionReturn(0);
194 }
195 EXTERN_C_END
196 
197 #undef __FUNCT__
198 #define __FUNCT__ "MatDestroy_MUMPS"
199 PetscErrorCode MatDestroy_MUMPS(Mat A)
200 {
201   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
202   PetscErrorCode ierr;
203   PetscMPIInt    size=lu->size;
204   PetscErrorCode (*specialdestroy)(Mat);
205   PetscFunctionBegin;
206   if (lu->CleanUpMUMPS) {
207     /* Terminate instance, deallocate memories */
208     if (size > 1){
209       ierr = PetscFree(lu->id.sol_loc);CHKERRQ(ierr);
210       ierr = VecScatterDestroy(lu->scat_rhs);CHKERRQ(ierr);
211       ierr = VecDestroy(lu->b_seq);CHKERRQ(ierr);
212       ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);
213       ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);
214       ierr = PetscFree(lu->val);CHKERRQ(ierr);
215     }
216     lu->id.job=JOB_END;
217 #if defined(PETSC_USE_COMPLEX)
218     zmumps_c(&lu->id);
219 #else
220     dmumps_c(&lu->id);
221 #endif
222     ierr = PetscFree(lu->irn);CHKERRQ(ierr);
223     ierr = PetscFree(lu->jcn);CHKERRQ(ierr);
224     ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr);
225   }
226   specialdestroy = lu->specialdestroy;
227   ierr = (*specialdestroy)(A);CHKERRQ(ierr);
228   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
229   PetscFunctionReturn(0);
230 }
231 
232 #undef __FUNCT__
233 #define __FUNCT__ "MatDestroy_AIJMUMPS"
234 PetscErrorCode MatDestroy_AIJMUMPS(Mat A)
235 {
236   PetscErrorCode ierr;
237   PetscMPIInt    size;
238 
239   PetscFunctionBegin;
240   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
241   if (size==1) {
242     ierr = MatConvert_MUMPS_Base(A,MATSEQAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
243   } else {
244     ierr = MatConvert_MUMPS_Base(A,MATMPIAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
245   }
246   PetscFunctionReturn(0);
247 }
248 
249 #undef __FUNCT__
250 #define __FUNCT__ "MatDestroy_SBAIJMUMPS"
251 PetscErrorCode MatDestroy_SBAIJMUMPS(Mat A)
252 {
253   PetscErrorCode ierr;
254   PetscMPIInt    size;
255 
256   PetscFunctionBegin;
257   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
258   if (size==1) {
259     ierr = MatConvert_MUMPS_Base(A,MATSEQSBAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
260   } else {
261     ierr = MatConvert_MUMPS_Base(A,MATMPISBAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
262   }
263   PetscFunctionReturn(0);
264 }
265 
266 #undef __FUNCT__
267 #define __FUNCT__ "MatSolve_MUMPS"
268 PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x) {
269   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
270   PetscScalar    *array;
271   Vec            x_seq;
272   IS             is_iden,is_petsc;
273   VecScatter     scat_rhs=lu->scat_rhs;
274   PetscErrorCode ierr;
275   PetscInt       i;
276 
277   PetscFunctionBegin;
278   lu->id.nrhs = 1;
279   x_seq = lu->b_seq;
280   if (lu->size > 1){
281     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
282     ierr = VecScatterBegin(scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
283     ierr = VecScatterEnd(scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
284     if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);}
285   } else {  /* size == 1 */
286     ierr = VecCopy(b,x);CHKERRQ(ierr);
287     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
288   }
289   if (!lu->myid) { /* define rhs on the host */
290 #if defined(PETSC_USE_COMPLEX)
291     lu->id.rhs = (mumps_double_complex*)array;
292 #else
293     lu->id.rhs = array;
294 #endif
295   }
296   if (lu->size == 1){
297     ierr = VecRestoreArray(x,&array);CHKERRQ(ierr);
298   } else if (!lu->myid){
299     ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr);
300   }
301 
302   if (lu->size > 1){
303     /* distributed solution */
304     lu->id.ICNTL(21) = 1;
305     if (!lu->nSolve){
306       /* Create x_seq=sol_loc for repeated use */
307       PetscInt    lsol_loc;
308       PetscScalar *sol_loc;
309       lsol_loc = lu->id.INFO(23); /* length of sol_loc */
310       ierr = PetscMalloc((1+lsol_loc)*(sizeof(PetscScalar)+sizeof(PetscInt)),&sol_loc);CHKERRQ(ierr);
311       lu->id.isol_loc = (PetscInt *)(sol_loc + lsol_loc);
312       lu->id.lsol_loc = lsol_loc;
313 #if defined(PETSC_USE_COMPLEX)
314       lu->id.sol_loc  = (ZMUMPS_DOUBLE *)sol_loc;
315 #else
316       lu->id.sol_loc  = (DMUMPS_DOUBLE *)sol_loc;
317 #endif
318       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,lsol_loc,sol_loc,&lu->x_seq);CHKERRQ(ierr);
319     }
320   }
321 
322   /* solve phase */
323   /*-------------*/
324   lu->id.job = 3;
325 #if defined(PETSC_USE_COMPLEX)
326   zmumps_c(&lu->id);
327 #else
328   dmumps_c(&lu->id);
329 #endif
330   if (lu->id.INFOG(1) < 0) {
331     SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1));
332   }
333 
334   if (lu->size > 1) { /* convert mumps distributed solution to petsc mpi x */
335     if (!lu->nSolve){ /* create scatter scat_sol */
336       ierr = ISCreateStride(PETSC_COMM_SELF,lu->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */
337       for (i=0; i<lu->id.lsol_loc; i++){
338         lu->id.isol_loc[i] -= 1; /* change Fortran style to C style */
339       }
340       ierr = ISCreateGeneral(PETSC_COMM_SELF,lu->id.lsol_loc,lu->id.isol_loc,&is_petsc);CHKERRQ(ierr);  /* to */
341       ierr = VecScatterCreate(lu->x_seq,is_iden,x,is_petsc,&lu->scat_sol);CHKERRQ(ierr);
342       ierr = ISDestroy(is_iden);CHKERRQ(ierr);
343       ierr = ISDestroy(is_petsc);CHKERRQ(ierr);
344     }
345     ierr = VecScatterBegin(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
346     ierr = VecScatterEnd(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
347   }
348   lu->nSolve++;
349   PetscFunctionReturn(0);
350 }
351 
352 /*
353   input:
354    F:        numeric factor
355   output:
356    nneg:     total number of negative pivots
357    nzero:    0
358    npos:     (global dimension of F) - nneg
359 */
360 
361 #undef __FUNCT__
362 #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
363 PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
364 {
365   Mat_MUMPS      *lu =(Mat_MUMPS*)F->spptr;
366   PetscErrorCode ierr;
367   PetscMPIInt    size;
368 
369   PetscFunctionBegin;
370   ierr = MPI_Comm_size(((PetscObject)F)->comm,&size);CHKERRQ(ierr);
371   /* 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 */
372   if (size > 1 && lu->id.ICNTL(13) != 1){
373     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));
374   }
375   if (nneg){
376     if (!lu->myid){
377       *nneg = lu->id.INFOG(12);
378     }
379     ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr);
380   }
381   if (nzero) *nzero = 0;
382   if (npos)  *npos  = F->rmap.N - (*nneg);
383   PetscFunctionReturn(0);
384 }
385 
386 #undef __FUNCT__
387 #define __FUNCT__ "MatFactorNumeric_MUMPS"
388 PetscErrorCode MatFactorNumeric_MUMPS(Mat A,MatFactorInfo *info,Mat *F)
389 {
390   Mat_MUMPS      *lu =(Mat_MUMPS*)(*F)->spptr;
391   Mat_MUMPS      *lua=(Mat_MUMPS*)(A)->spptr;
392   PetscErrorCode ierr;
393   PetscInt       rnz,nnz,nz=0,i,M=A->rmap.N,*ai,*aj,icntl;
394   PetscTruth     valOnly,flg;
395   Mat            F_diag;
396 
397   PetscFunctionBegin;
398   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
399     (*F)->ops->solve    = MatSolve_MUMPS;
400 
401     /* Initialize a MUMPS instance */
402     ierr = MPI_Comm_rank(((PetscObject)A)->comm, &lu->myid);
403     ierr = MPI_Comm_size(((PetscObject)A)->comm,&lu->size);CHKERRQ(ierr);
404     lua->myid = lu->myid; lua->size = lu->size;
405     lu->id.job = JOB_INIT;
406     ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->comm_mumps));CHKERRQ(ierr);
407     ierr = MPICCommToFortranComm(lu->comm_mumps,&(lu->id.comm_fortran));CHKERRQ(ierr);
408 
409     /* Set mumps options */
410     ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
411     lu->id.par=1;  /* host participates factorizaton and solve */
412     lu->id.sym=lu->sym;
413     if (lu->sym == 2){
414       ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr);
415       if (flg && icntl == 1) lu->id.sym=icntl;  /* matrix is spd */
416     }
417 #if defined(PETSC_USE_COMPLEX)
418     zmumps_c(&lu->id);
419 #else
420     dmumps_c(&lu->id);
421 #endif
422 
423     if (lu->size == 1){
424       lu->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
425     } else {
426       lu->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
427     }
428 
429     icntl=-1;
430     lu->id.ICNTL(4) = 0;  /* level of printing; overwrite mumps default ICNTL(4)=2 */
431     ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
432     if ((flg && icntl > 0) || PetscLogPrintInfo) {
433       lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */
434     } else { /* no output */
435       lu->id.ICNTL(1) = 0;  /* error message, default= 6 */
436       lu->id.ICNTL(2) = -1; /* output stream for diagnostic printing, statistics, and warning. default=0 */
437       lu->id.ICNTL(3) = -1; /* output stream for global information, default=6 */
438     }
439     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);
440     icntl=-1;
441     ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
442     if (flg) {
443       if (icntl== 1){
444         SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n");
445       } else {
446         lu->id.ICNTL(7) = icntl;
447       }
448     }
449     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);
450     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);
451     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);
452     ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr);
453     ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr);
454     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);
455     ierr = PetscOptionsInt("-mat_mumps_icntl_15","ICNTL(15): efficiency control","None",lu->id.ICNTL(15),&lu->id.ICNTL(15),PETSC_NULL);CHKERRQ(ierr);
456 
457     ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",lu->id.CNTL(1),&lu->id.CNTL(1),PETSC_NULL);CHKERRQ(ierr);
458     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);
459     ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",lu->id.CNTL(3),&lu->id.CNTL(3),PETSC_NULL);CHKERRQ(ierr);
460     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);
461     PetscOptionsEnd();
462   }
463 
464   /* define matrix A */
465   switch (lu->id.ICNTL(18)){
466   case 0:  /* centralized assembled matrix input (size=1) */
467     if (!lu->myid) {
468       if (lua->isAIJ){
469         Mat_SeqAIJ   *aa = (Mat_SeqAIJ*)A->data;
470         nz               = aa->nz;
471         ai = aa->i; aj = aa->j; lu->val = aa->a;
472       } else {
473         Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data;
474         nz                  =  aa->nz;
475         ai = aa->i; aj = aa->j; lu->val = aa->a;
476       }
477       if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */
478         ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->irn);CHKERRQ(ierr);
479         ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->jcn);CHKERRQ(ierr);
480         nz = 0;
481         for (i=0; i<M; i++){
482           rnz = ai[i+1] - ai[i];
483           while (rnz--) {  /* Fortran row/col index! */
484             lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++;
485           }
486         }
487       }
488     }
489     break;
490   case 3:  /* distributed assembled matrix input (size>1) */
491     if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
492       valOnly = PETSC_FALSE;
493     } else {
494       valOnly = PETSC_TRUE; /* only update mat values, not row and col index */
495     }
496     ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
497     break;
498   default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS.");
499   }
500 
501   /* analysis phase */
502   /*----------------*/
503   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
504     lu->id.job = 1;
505 
506     lu->id.n = M;
507     switch (lu->id.ICNTL(18)){
508     case 0:  /* centralized assembled matrix input */
509       if (!lu->myid) {
510         lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
511         if (lu->id.ICNTL(6)>1){
512 #if defined(PETSC_USE_COMPLEX)
513           lu->id.a = (mumps_double_complex*)lu->val;
514 #else
515           lu->id.a = lu->val;
516 #endif
517         }
518       }
519       break;
520     case 3:  /* distributed assembled matrix input (size>1) */
521       lu->id.nz_loc = nnz;
522       lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
523       if (lu->id.ICNTL(6)>1) {
524 #if defined(PETSC_USE_COMPLEX)
525         lu->id.a_loc = (mumps_double_complex*)lu->val;
526 #else
527         lu->id.a_loc = lu->val;
528 #endif
529       }
530       /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
531       IS  is_iden;
532       Vec b;
533       if (!lu->myid){
534         ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap.N,&lu->b_seq);CHKERRQ(ierr);
535         ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap.N,0,1,&is_iden);CHKERRQ(ierr);
536       } else {
537         ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr);
538         ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
539       }
540       ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr);
541       ierr = VecSetSizes(b,A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr);
542       ierr = VecSetFromOptions(b);CHKERRQ(ierr);
543 
544       ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr);
545       ierr = ISDestroy(is_iden);CHKERRQ(ierr);
546       ierr = VecDestroy(b);CHKERRQ(ierr);
547       break;
548     }
549 #if defined(PETSC_USE_COMPLEX)
550     zmumps_c(&lu->id);
551 #else
552     dmumps_c(&lu->id);
553 #endif
554     if (lu->id.INFOG(1) < 0) {
555       SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
556     }
557   }
558 
559   /* numerical factorization phase */
560   /*-------------------------------*/
561   lu->id.job = 2;
562   if(!lu->id.ICNTL(18)) {
563     if (!lu->myid) {
564 #if defined(PETSC_USE_COMPLEX)
565       lu->id.a = (mumps_double_complex*)lu->val;
566 #else
567       lu->id.a = lu->val;
568 #endif
569     }
570   } else {
571 #if defined(PETSC_USE_COMPLEX)
572     lu->id.a_loc = (mumps_double_complex*)lu->val;
573 #else
574     lu->id.a_loc = lu->val;
575 #endif
576   }
577 #if defined(PETSC_USE_COMPLEX)
578   zmumps_c(&lu->id);
579 #else
580   dmumps_c(&lu->id);
581 #endif
582   if (lu->id.INFOG(1) < 0) {
583     if (lu->id.INFO(1) == -13) {
584       SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",lu->id.INFO(2));
585     } else {
586       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));
587     }
588   }
589 
590   if (!lu->myid && lu->id.ICNTL(16) > 0){
591     SETERRQ1(PETSC_ERR_LIB,"  lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16));
592   }
593 
594   if (lu->size > 1){
595     if ((*F)->factor == FACTOR_LU){
596       F_diag = ((Mat_MPIAIJ *)(*F)->data)->A;
597     } else {
598       F_diag = ((Mat_MPISBAIJ *)(*F)->data)->A;
599     }
600     F_diag->assembled = PETSC_TRUE;
601     if (lu->nSolve){
602       ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);
603       ierr = PetscFree(lu->id.sol_loc);CHKERRQ(ierr);
604       ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);
605     }
606   }
607   (*F)->assembled   = PETSC_TRUE;
608   lu->matstruc      = SAME_NONZERO_PATTERN;
609   lu->CleanUpMUMPS  = PETSC_TRUE;
610   lu->nSolve        = 0;
611   PetscFunctionReturn(0);
612 }
613 
614 /* Note the Petsc r and c permutations are ignored */
615 #undef __FUNCT__
616 #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
617 PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) {
618   Mat            B;
619   Mat_MUMPS      *lu;
620   PetscErrorCode ierr;
621 
622   PetscFunctionBegin;
623   /* Create the factorization matrix */
624   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
625   ierr = MatSetSizes(B,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);CHKERRQ(ierr);
626   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
627   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
628   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
629 
630   B->ops->lufactornumeric = MatFactorNumeric_MUMPS;
631   B->factor               = FACTOR_LU;
632   lu                      = (Mat_MUMPS*)B->spptr;
633   lu->sym                 = 0;
634   lu->matstruc            = DIFFERENT_NONZERO_PATTERN;
635 
636   *F = B;
637   PetscFunctionReturn(0);
638 }
639 
640 /* Note the Petsc r permutation is ignored */
641 #undef __FUNCT__
642 #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS"
643 PetscErrorCode MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat A,IS r,MatFactorInfo *info,Mat *F) {
644   Mat            B;
645   Mat_MUMPS      *lu;
646   PetscErrorCode ierr;
647 
648   PetscFunctionBegin;
649   /* Create the factorization matrix */
650   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
651   ierr = MatSetSizes(B,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);CHKERRQ(ierr);
652   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
653   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
654   ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
655 
656   B->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
657   B->ops->getinertia            = MatGetInertia_SBAIJMUMPS;
658   B->factor                     = FACTOR_CHOLESKY;
659   lu                            = (Mat_MUMPS*)B->spptr;
660   lu->sym                       = 2;
661   lu->matstruc                  = DIFFERENT_NONZERO_PATTERN;
662 
663   *F = B;
664   PetscFunctionReturn(0);
665 }
666 
667 #undef __FUNCT__
668 #define __FUNCT__ "MatFactorInfo_MUMPS"
669 PetscErrorCode MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) {
670   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
671   PetscErrorCode ierr;
672 
673   PetscFunctionBegin;
674   /* check if matrix is mumps type */
675   if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0);
676 
677   ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
678   ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                  %d \n",lu->id.sym);CHKERRQ(ierr);
679   ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):           %d \n",lu->id.par);CHKERRQ(ierr);
680   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):        %d \n",lu->id.ICNTL(1));CHKERRQ(ierr);
681   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg):%d \n",lu->id.ICNTL(2));CHKERRQ(ierr);
682   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):  %d \n",lu->id.ICNTL(3));CHKERRQ(ierr);
683   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):       %d \n",lu->id.ICNTL(4));CHKERRQ(ierr);
684   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):        %d \n",lu->id.ICNTL(5));CHKERRQ(ierr);
685   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):       %d \n",lu->id.ICNTL(6));CHKERRQ(ierr);
686   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (matrix ordering):         %d \n",lu->id.ICNTL(7));CHKERRQ(ierr);
687   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scalling strategy):       %d \n",lu->id.ICNTL(8));CHKERRQ(ierr);
688   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(9) (A/A^T x=b is solved):     %d \n",lu->id.ICNTL(9));CHKERRQ(ierr);
689   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr);
690   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):         %d \n",lu->id.ICNTL(11));CHKERRQ(ierr);
691   if (!lu->myid && lu->id.ICNTL(11)>0) {
692     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(4) (inf norm of input mat):        %g\n",lu->id.RINFOG(4));CHKERRQ(ierr);
693     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(5) (inf norm of solution):         %g\n",lu->id.RINFOG(5));CHKERRQ(ierr);
694     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(6) (inf norm of residual):         %g\n",lu->id.RINFOG(6));CHKERRQ(ierr);
695     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);
696     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(9) (error estimate):               %g \n",lu->id.RINFOG(9));CHKERRQ(ierr);
697     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",lu->id.RINFOG(10),lu->id.RINFOG(11));CHKERRQ(ierr);
698 
699   }
700   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",lu->id.ICNTL(12));CHKERRQ(ierr);
701   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",lu->id.ICNTL(13));CHKERRQ(ierr);
702   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr);
703   /* ICNTL(15-17) not used */
704   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",lu->id.ICNTL(18));CHKERRQ(ierr);
705   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Shur complement info):                       %d \n",lu->id.ICNTL(19));CHKERRQ(ierr);
706   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",lu->id.ICNTL(20));CHKERRQ(ierr);
707   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d \n",lu->id.ICNTL(21));CHKERRQ(ierr);
708 
709   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",lu->id.CNTL(1));CHKERRQ(ierr);
710   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr);
711   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",lu->id.CNTL(3));CHKERRQ(ierr);
712   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",lu->id.CNTL(4));CHKERRQ(ierr);
713 
714   /* infomation local to each processor */
715   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);}
716   ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm,"             [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr);
717   ierr = PetscSynchronizedFlush(((PetscObject)A)->comm);
718   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);}
719   ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr);
720   ierr = PetscSynchronizedFlush(((PetscObject)A)->comm);
721   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);}
722   ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr);
723   ierr = PetscSynchronizedFlush(((PetscObject)A)->comm);
724   /*
725   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      INFO(2) (info about error or warning ): \n");CHKERRQ(ierr);}
726   ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(2));CHKERRQ(ierr);
727   ierr = PetscSynchronizedFlush(((PetscObject)A)->comm);
728   */
729 
730   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr);}
731   ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(15));CHKERRQ(ierr);
732   ierr = PetscSynchronizedFlush(((PetscObject)A)->comm);
733 
734   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr);}
735   ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(16));CHKERRQ(ierr);
736   ierr = PetscSynchronizedFlush(((PetscObject)A)->comm);
737 
738   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);}
739   ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(23));CHKERRQ(ierr);
740   ierr = PetscSynchronizedFlush(((PetscObject)A)->comm);
741 
742   if (!lu->myid){ /* information from the host */
743     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr);
744     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr);
745     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr);
746 
747     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr);
748     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr);
749     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr);
750     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr);
751     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr);
752     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr);
753     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",lu->id.INFOG(9));CHKERRQ(ierr);
754     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",lu->id.INFOG(10));CHKERRQ(ierr);
755     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",lu->id.INFOG(11));CHKERRQ(ierr);
756     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr);
757     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr);
758     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr);
759     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr);
760     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);
761     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);
762     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);
763     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);
764      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr);
765      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);
766      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);
767      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",lu->id.INFOG(23));CHKERRQ(ierr);
768      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",lu->id.INFOG(24));CHKERRQ(ierr);
769      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",lu->id.INFOG(25));CHKERRQ(ierr);
770   }
771 
772   PetscFunctionReturn(0);
773 }
774 
775 #undef __FUNCT__
776 #define __FUNCT__ "MatView_MUMPS"
777 PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer) {
778   PetscErrorCode    ierr;
779   PetscTruth        iascii;
780   PetscViewerFormat format;
781   Mat_MUMPS         *mumps=(Mat_MUMPS*)(A->spptr);
782 
783   PetscFunctionBegin;
784   ierr = (*mumps->MatView)(A,viewer);CHKERRQ(ierr);
785 
786   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
787   if (iascii) {
788     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
789     if (format == PETSC_VIEWER_ASCII_INFO){
790       ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr);
791     }
792   }
793   PetscFunctionReturn(0);
794 }
795 
796 #undef __FUNCT__
797 #define __FUNCT__ "MatAssemblyEnd_AIJMUMPS"
798 PetscErrorCode MatAssemblyEnd_AIJMUMPS(Mat A,MatAssemblyType mode) {
799   PetscErrorCode ierr;
800   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
801 
802   PetscFunctionBegin;
803   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
804 
805   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
806   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
807   A->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
808   PetscFunctionReturn(0);
809 }
810 
811 EXTERN_C_BEGIN
812 #undef __FUNCT__
813 #define __FUNCT__ "MatConvert_AIJ_AIJMUMPS"
814 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_AIJ_AIJMUMPS(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
815 {
816   PetscErrorCode ierr;
817   PetscMPIInt    size;
818   MPI_Comm       comm;
819   Mat            B=*newmat;
820   Mat_MUMPS      *mumps;
821 
822   PetscFunctionBegin;
823   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
824   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
825 
826   if (reuse == MAT_INITIAL_MATRIX) {
827     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
828     /* A may have special container that is not duplicated,
829        e.g., A is obtainted from MatMatMult(,&A). Save B->ops instead */
830     mumps->MatDuplicate              = B->ops->duplicate;
831     mumps->MatView                   = B->ops->view;
832     mumps->MatAssemblyEnd            = B->ops->assemblyend;
833     mumps->MatLUFactorSymbolic       = B->ops->lufactorsymbolic;
834     mumps->MatCholeskyFactorSymbolic = B->ops->choleskyfactorsymbolic;
835     mumps->MatDestroy                = B->ops->destroy;
836   } else {
837     mumps->MatDuplicate              = A->ops->duplicate;
838     mumps->MatView                   = A->ops->view;
839     mumps->MatAssemblyEnd            = A->ops->assemblyend;
840     mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
841     mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
842     mumps->MatDestroy                = A->ops->destroy;
843   }
844   mumps->specialdestroy            = MatDestroy_AIJMUMPS;
845   mumps->CleanUpMUMPS              = PETSC_FALSE;
846   mumps->isAIJ                     = PETSC_TRUE;
847 
848   B->spptr                         = (void*)mumps;
849   B->ops->duplicate                = MatDuplicate_MUMPS;
850   B->ops->view                     = MatView_MUMPS;
851   B->ops->assemblyend              = MatAssemblyEnd_AIJMUMPS;
852   B->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
853   B->ops->destroy                  = MatDestroy_MUMPS;
854 
855   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
856   if (size == 1) {
857     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_aijmumps_C",
858                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
859     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_seqaij_C",
860                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
861   } else {
862     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpiaij_aijmumps_C",
863                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
864     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_mpiaij_C",
865                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
866   }
867 
868   ierr = PetscInfo(0,"Using MUMPS for LU factorization and solves.\n");CHKERRQ(ierr);
869   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
870   *newmat = B;
871   PetscFunctionReturn(0);
872 }
873 EXTERN_C_END
874 
875 /*MC
876   MATAIJMUMPS - MATAIJMUMPS = "aijmumps" - A matrix type providing direct solvers (LU) for distributed
877   and sequential matrices via the external package MUMPS.
878 
879   If MUMPS is installed (see the manual for instructions
880   on how to declare the existence of external packages),
881   a matrix type can be constructed which invokes MUMPS solvers.
882   After calling MatCreate(...,A), simply call MatSetType(A,MATAIJMUMPS), then
883   optionally call MatSeqAIJSetPreallocation() or MatMPIAIJSetPreallocation() etc DO NOT
884   call MatCreateSeqAIJ/MPIAIJ() directly or the preallocation information will be LOST!
885 
886   If created with a single process communicator, this matrix type inherits from MATSEQAIJ.
887   Otherwise, this matrix type inherits from MATMPIAIJ.  Hence for single process communicators,
888   MatSeqAIJSetPreallocation() is supported, and similarly MatMPIAIJSetPreallocation() is supported
889   for communicators controlling multiple processes.  It is recommended that you call both of
890   the above preallocation routines for simplicity.  One can also call MatConvert() for an inplace
891   conversion to or from the MATSEQAIJ or MATMPIAIJ type (depending on the communicator size)
892   without data copy AFTER the matrix values are set.
893 
894   Options Database Keys:
895 + -mat_type aijmumps - sets the matrix type to "aijmumps" during a call to MatSetFromOptions()
896 . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
897 . -mat_mumps_icntl_4 <0,1,2,3,4> - print level
898 . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
899 . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
900 . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
901 . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
902 . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
903 . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
904 . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
905 . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
906 . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
907 . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
908 . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
909 - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
910 
911   Level: beginner
912 
913 .seealso: MATSBAIJMUMPS
914 M*/
915 
916 EXTERN_C_BEGIN
917 #undef __FUNCT__
918 #define __FUNCT__ "MatCreate_AIJMUMPS"
919 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_AIJMUMPS(Mat A)
920 {
921   PetscErrorCode ierr;
922   PetscMPIInt    size;
923 
924   PetscFunctionBegin;
925   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
926   if (size == 1) {
927     ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
928   } else {
929     ierr   = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
930     /*
931     Mat A_diag = ((Mat_MPIAIJ *)A->data)->A;
932     ierr   = MatConvert_AIJ_AIJMUMPS(A_diag,MATAIJMUMPS,MAT_REUSE_MATRIX,&A_diag);CHKERRQ(ierr);
933     */
934   }
935   ierr = MatConvert_AIJ_AIJMUMPS(A,MATAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
936   PetscFunctionReturn(0);
937 }
938 EXTERN_C_END
939 
940 #undef __FUNCT__
941 #define __FUNCT__ "MatAssemblyEnd_SBAIJMUMPS"
942 PetscErrorCode MatAssemblyEnd_SBAIJMUMPS(Mat A,MatAssemblyType mode)
943 {
944   PetscErrorCode ierr;
945   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
946 
947   PetscFunctionBegin;
948   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
949   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
950   A->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
951   PetscFunctionReturn(0);
952 }
953 
954 EXTERN_C_BEGIN
955 #undef __FUNCT__
956 #define __FUNCT__ "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS"
957 PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation_MPISBAIJMUMPS(Mat B,PetscInt bs,PetscInt d_nz,PetscInt *d_nnz,PetscInt o_nz,PetscInt *o_nnz)
958 {
959   Mat       A;
960   Mat_MUMPS *mumps=(Mat_MUMPS*)B->spptr;
961   PetscErrorCode ierr;
962 
963   PetscFunctionBegin;
964   /*
965     After performing the MPISBAIJ Preallocation, we need to convert the local diagonal block matrix
966     into MUMPS type so that the block jacobi preconditioner (for example) can use MUMPS.  I would
967     like this to be done in the MatCreate routine, but the creation of this inner matrix requires
968     block size info so that PETSc can determine the local size properly.  The block size info is set
969     in the preallocation routine.
970   */
971   ierr = (*mumps->MatPreallocate)(B,bs,d_nz,d_nnz,o_nz,o_nnz);
972   A    = ((Mat_MPISBAIJ *)B->data)->A;
973   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
974   PetscFunctionReturn(0);
975 }
976 EXTERN_C_END
977 
978 EXTERN_C_BEGIN
979 #undef __FUNCT__
980 #define __FUNCT__ "MatConvert_SBAIJ_SBAIJMUMPS"
981 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SBAIJ_SBAIJMUMPS(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
982 {
983   PetscErrorCode ierr;
984   PetscMPIInt    size;
985   MPI_Comm       comm;
986   Mat            B=*newmat;
987   Mat_MUMPS      *mumps;
988   void           (*f)(void);
989 
990   PetscFunctionBegin;
991   if (reuse == MAT_INITIAL_MATRIX) {
992     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
993   }
994 
995   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
996   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
997 
998   mumps->MatDuplicate              = A->ops->duplicate;
999   mumps->MatView                   = A->ops->view;
1000   mumps->MatAssemblyEnd            = A->ops->assemblyend;
1001   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
1002   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
1003   mumps->MatDestroy                = A->ops->destroy;
1004   mumps->specialdestroy            = MatDestroy_SBAIJMUMPS;
1005   mumps->CleanUpMUMPS              = PETSC_FALSE;
1006   mumps->isAIJ                     = PETSC_FALSE;
1007 
1008   B->spptr                         = (void*)mumps;
1009   B->ops->duplicate                = MatDuplicate_MUMPS;
1010   B->ops->view                     = MatView_MUMPS;
1011   B->ops->assemblyend              = MatAssemblyEnd_SBAIJMUMPS;
1012   B->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
1013   B->ops->destroy                  = MatDestroy_MUMPS;
1014 
1015   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
1016   if (size == 1) {
1017     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C",
1018                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
1019     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C",
1020                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
1021   } else {
1022   /* I really don't like needing to know the tag: MatMPISBAIJSetPreallocation_C */
1023     ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr);
1024     if (f) { /* This case should always be true when this routine is called */
1025       mumps->MatPreallocate = (PetscErrorCode (*)(Mat,int,int,int*,int,int*))f;
1026       ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPISBAIJSetPreallocation_C",
1027                                                "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS",
1028                                                MatMPISBAIJSetPreallocation_MPISBAIJMUMPS);CHKERRQ(ierr);
1029     }
1030     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C",
1031                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
1032     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C",
1033                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
1034   }
1035 
1036   ierr = PetscInfo(0,"Using MUMPS for Cholesky factorization and solves.\n");CHKERRQ(ierr);
1037   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
1038   *newmat = B;
1039   PetscFunctionReturn(0);
1040 }
1041 EXTERN_C_END
1042 
1043 #undef __FUNCT__
1044 #define __FUNCT__ "MatDuplicate_MUMPS"
1045 PetscErrorCode MatDuplicate_MUMPS(Mat A, MatDuplicateOption op, Mat *M) {
1046   PetscErrorCode ierr;
1047   Mat_MUMPS   *lu=(Mat_MUMPS *)A->spptr;
1048 
1049   PetscFunctionBegin;
1050   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
1051   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_MUMPS));CHKERRQ(ierr);
1052   PetscFunctionReturn(0);
1053 }
1054 
1055 /*MC
1056   MATSBAIJMUMPS - MATSBAIJMUMPS = "sbaijmumps" - A symmetric matrix type providing direct solvers (Cholesky) for
1057   distributed and sequential matrices via the external package MUMPS.
1058 
1059   If MUMPS is installed (see the manual for instructions
1060   on how to declare the existence of external packages),
1061   a matrix type can be constructed which invokes MUMPS solvers.
1062   After calling MatCreate(...,A), simply call MatSetType(A,MATSBAIJMUMPS), then
1063   optionally call MatSeqSBAIJSetPreallocation() or MatMPISBAIJSetPreallocation() DO NOT
1064   call MatCreateSeqSBAIJ/MPISBAIJ() directly or the preallocation information will be LOST!
1065 
1066   If created with a single process communicator, this matrix type inherits from MATSEQSBAIJ.
1067   Otherwise, this matrix type inherits from MATMPISBAIJ.  Hence for single process communicators,
1068   MatSeqSBAIJSetPreallocation() is supported, and similarly MatMPISBAIJSetPreallocation() is supported
1069   for communicators controlling multiple processes.  It is recommended that you call both of
1070   the above preallocation routines for simplicity.  One can also call MatConvert() for an inplace
1071   conversion to or from the MATSEQSBAIJ or MATMPISBAIJ type (depending on the communicator size)
1072   without data copy AFTER the matrix values have been set.
1073 
1074   Options Database Keys:
1075 + -mat_type sbaijmumps - sets the matrix type to "sbaijmumps" during a call to MatSetFromOptions()
1076 . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
1077 . -mat_mumps_icntl_4 <0,...,4> - print level
1078 . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
1079 . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
1080 . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
1081 . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
1082 . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
1083 . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
1084 . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
1085 . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
1086 . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
1087 . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
1088 . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
1089 - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
1090 
1091   Level: beginner
1092 
1093 .seealso: MATAIJMUMPS
1094 M*/
1095 
1096 EXTERN_C_BEGIN
1097 #undef __FUNCT__
1098 #define __FUNCT__ "MatCreate_SBAIJMUMPS"
1099 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SBAIJMUMPS(Mat A)
1100 {
1101   PetscErrorCode ierr;
1102   PetscMPIInt    size;
1103 
1104   PetscFunctionBegin;
1105   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
1106   if (size == 1) {
1107     ierr = MatSetType(A,MATSEQSBAIJ);CHKERRQ(ierr);
1108   } else {
1109     ierr   = MatSetType(A,MATMPISBAIJ);CHKERRQ(ierr);
1110   }
1111   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
1112   PetscFunctionReturn(0);
1113 }
1114 EXTERN_C_END
1115