xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision 9481e6e9e9fee918c57079e292375d33e010da0c)
1 
2 /*
3     Provides an interface to the MUMPS sparse solver
4 */
5 
6 #include <../src/mat/impls/aij/mpi/mpiaij.h> /*I  "petscmat.h"  I*/
7 #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
8 
9 EXTERN_C_BEGIN
10 #if defined(PETSC_USE_COMPLEX)
11 #if defined(PETSC_USE_REAL_SINGLE)
12 #include <cmumps_c.h>
13 #else
14 #include <zmumps_c.h>
15 #endif
16 #else
17 #if defined(PETSC_USE_REAL_SINGLE)
18 #include <smumps_c.h>
19 #else
20 #include <dmumps_c.h>
21 #endif
22 #endif
23 EXTERN_C_END
24 #define JOB_INIT -1
25 #define JOB_FACTSYMBOLIC 1
26 #define JOB_FACTNUMERIC 2
27 #define JOB_SOLVE 3
28 #define JOB_END -2
29 
30 /* calls to MUMPS */
31 #if defined(PETSC_USE_COMPLEX)
32 #if defined(PETSC_USE_REAL_SINGLE)
33 #define PetscMUMPS_c cmumps_c
34 #else
35 #define PetscMUMPS_c zmumps_c
36 #endif
37 #else
38 #if defined(PETSC_USE_REAL_SINGLE)
39 #define PetscMUMPS_c smumps_c
40 #else
41 #define PetscMUMPS_c dmumps_c
42 #endif
43 #endif
44 
45 /* declare MumpsScalar */
46 #if defined(PETSC_USE_COMPLEX)
47 #if defined(PETSC_USE_REAL_SINGLE)
48 #define MumpsScalar mumps_complex
49 #else
50 #define MumpsScalar mumps_double_complex
51 #endif
52 #else
53 #define MumpsScalar PetscScalar
54 #endif
55 
56 /* macros s.t. indices match MUMPS documentation */
57 #define ICNTL(I) icntl[(I)-1]
58 #define CNTL(I) cntl[(I)-1]
59 #define INFOG(I) infog[(I)-1]
60 #define INFO(I) info[(I)-1]
61 #define RINFOG(I) rinfog[(I)-1]
62 #define RINFO(I) rinfo[(I)-1]
63 
64 typedef struct {
65 #if defined(PETSC_USE_COMPLEX)
66 #if defined(PETSC_USE_REAL_SINGLE)
67   CMUMPS_STRUC_C id;
68 #else
69   ZMUMPS_STRUC_C id;
70 #endif
71 #else
72 #if defined(PETSC_USE_REAL_SINGLE)
73   SMUMPS_STRUC_C id;
74 #else
75   DMUMPS_STRUC_C id;
76 #endif
77 #endif
78 
79   MatStructure matstruc;
80   PetscMPIInt  myid,size;
81   PetscInt     *irn,*jcn,nz,sym;
82   PetscScalar  *val;
83   MPI_Comm     comm_mumps;
84   PetscBool    isAIJ;
85   PetscInt     ICNTL9_pre;           /* check if ICNTL(9) is changed from previous MatSolve */
86   VecScatter   scat_rhs, scat_sol;   /* used by MatSolve() */
87   Vec          b_seq,x_seq;
88   PetscInt     ninfo,*info;          /* display INFO */
89   PetscInt     sizeredrhs;
90   PetscScalar  *schur_sol;
91   PetscInt     schur_sizesol;
92 
93   PetscErrorCode (*ConvertToTriples)(Mat, int, MatReuse, int*, int**, int**, PetscScalar**);
94 } Mat_MUMPS;
95 
96 extern PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);
97 
98 static PetscErrorCode MatMumpsResetSchur_Private(Mat_MUMPS* mumps)
99 {
100   PetscErrorCode ierr;
101 
102   PetscFunctionBegin;
103   ierr = PetscFree2(mumps->id.listvar_schur,mumps->id.schur);CHKERRQ(ierr);
104   ierr = PetscFree(mumps->id.redrhs);CHKERRQ(ierr);
105   ierr = PetscFree(mumps->schur_sol);CHKERRQ(ierr);
106   mumps->id.size_schur = 0;
107   mumps->id.schur_lld  = 0;
108   mumps->id.ICNTL(19)  = 0;
109   PetscFunctionReturn(0);
110 }
111 
112 /* solve with rhs in mumps->id.redrhs and return in the same location */
113 static PetscErrorCode MatMumpsSolveSchur_Private(Mat F)
114 {
115   Mat_MUMPS            *mumps=(Mat_MUMPS*)F->data;
116   Mat                  S,B,X;
117   MatFactorSchurStatus schurstatus;
118   PetscInt             sizesol;
119   PetscErrorCode       ierr;
120 
121   PetscFunctionBegin;
122   ierr = MatFactorFactorizeSchurComplement(F);CHKERRQ(ierr);
123   ierr = MatFactorGetSchurComplement(F,&S,&schurstatus);CHKERRQ(ierr);
124   ierr = MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.nrhs,(PetscScalar*)mumps->id.redrhs,&B);CHKERRQ(ierr);
125   switch (schurstatus) {
126   case MAT_FACTOR_SCHUR_FACTORED:
127     ierr = MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.nrhs,(PetscScalar*)mumps->id.redrhs,&X);CHKERRQ(ierr);
128     if (!mumps->id.ICNTL(9)) { /* transpose solve */
129       ierr = MatMatSolveTranspose(S,B,X);CHKERRQ(ierr);
130     } else {
131       ierr = MatMatSolve(S,B,X);CHKERRQ(ierr);
132     }
133     break;
134   case MAT_FACTOR_SCHUR_INVERTED:
135     sizesol = mumps->id.nrhs*mumps->id.size_schur;
136     if (!mumps->schur_sol || sizesol > mumps->schur_sizesol) {
137       ierr = PetscFree(mumps->schur_sol);CHKERRQ(ierr);
138       ierr = PetscMalloc1(sizesol,&mumps->schur_sol);CHKERRQ(ierr);
139       mumps->schur_sizesol = sizesol;
140     }
141     ierr = MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.nrhs,mumps->schur_sol,&X);CHKERRQ(ierr);
142     if (!mumps->id.ICNTL(9)) { /* transpose solve */
143       ierr = MatTransposeMatMult(S,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);CHKERRQ(ierr);
144     } else {
145       ierr = MatMatMult(S,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);CHKERRQ(ierr);
146     }
147     ierr = MatCopy(X,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
148     break;
149   default:
150     SETERRQ1(PetscObjectComm((PetscObject)F),PETSC_ERR_SUP,"Unhandled MatFactorSchurStatus %D",F->schur_status);
151     break;
152   }
153   ierr = MatFactorRestoreSchurComplement(F,&S,schurstatus);CHKERRQ(ierr);
154   ierr = MatDestroy(&B);CHKERRQ(ierr);
155   ierr = MatDestroy(&X);CHKERRQ(ierr);
156   PetscFunctionReturn(0);
157 }
158 
159 static PetscErrorCode MatMumpsHandleSchur_Private(Mat F, PetscBool expansion)
160 {
161   Mat_MUMPS     *mumps=(Mat_MUMPS*)F->data;
162   PetscErrorCode ierr;
163 
164   PetscFunctionBegin;
165   if (!mumps->id.ICNTL(19)) { /* do nothing when Schur complement has not been computed */
166     PetscFunctionReturn(0);
167   }
168   if (!expansion) { /* prepare for the condensation step */
169     PetscInt sizeredrhs = mumps->id.nrhs*mumps->id.size_schur;
170     /* allocate MUMPS internal array to store reduced right-hand sides */
171     if (!mumps->id.redrhs || sizeredrhs > mumps->sizeredrhs) {
172       ierr = PetscFree(mumps->id.redrhs);CHKERRQ(ierr);
173       mumps->id.lredrhs = mumps->id.size_schur;
174       ierr = PetscMalloc1(mumps->id.nrhs*mumps->id.lredrhs,&mumps->id.redrhs);CHKERRQ(ierr);
175       mumps->sizeredrhs = mumps->id.nrhs*mumps->id.lredrhs;
176     }
177     mumps->id.ICNTL(26) = 1; /* condensation phase */
178   } else { /* prepare for the expansion step */
179     /* solve Schur complement (this has to be done by the MUMPS user, so basically us) */
180     ierr = MatMumpsSolveSchur_Private(F);CHKERRQ(ierr);
181     mumps->id.ICNTL(26) = 2; /* expansion phase */
182     PetscMUMPS_c(&mumps->id);
183     if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
184     /* restore defaults */
185     mumps->id.ICNTL(26) = -1;
186     /* free MUMPS internal array for redrhs if we have solved for multiple rhs in order to save memory space */
187     if (mumps->id.nrhs > 1) {
188       ierr = PetscFree(mumps->id.redrhs);CHKERRQ(ierr);
189       mumps->id.lredrhs = 0;
190       mumps->sizeredrhs = 0;
191     }
192   }
193   PetscFunctionReturn(0);
194 }
195 
196 /*
197   MatConvertToTriples_A_B - convert Petsc matrix to triples: row[nz], col[nz], val[nz]
198 
199   input:
200     A       - matrix in aij,baij or sbaij (bs=1) format
201     shift   - 0: C style output triple; 1: Fortran style output triple.
202     reuse   - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple
203               MAT_REUSE_MATRIX:   only the values in v array are updated
204   output:
205     nnz     - dim of r, c, and v (number of local nonzero entries of A)
206     r, c, v - row and col index, matrix values (matrix triples)
207 
208   The returned values r, c, and sometimes v are obtained in a single PetscMalloc(). Then in MatDestroy_MUMPS() it is
209   freed with PetscFree((mumps->irn);  This is not ideal code, the fact that v is ONLY sometimes part of mumps->irn means
210   that the PetscMalloc() cannot easily be replaced with a PetscMalloc3().
211 
212  */
213 
214 PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
215 {
216   const PetscInt *ai,*aj,*ajj,M=A->rmap->n;
217   PetscInt       nz,rnz,i,j;
218   PetscErrorCode ierr;
219   PetscInt       *row,*col;
220   Mat_SeqAIJ     *aa=(Mat_SeqAIJ*)A->data;
221 
222   PetscFunctionBegin;
223   *v=aa->a;
224   if (reuse == MAT_INITIAL_MATRIX) {
225     nz   = aa->nz;
226     ai   = aa->i;
227     aj   = aa->j;
228     *nnz = nz;
229     ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr);
230     col  = row + nz;
231 
232     nz = 0;
233     for (i=0; i<M; i++) {
234       rnz = ai[i+1] - ai[i];
235       ajj = aj + ai[i];
236       for (j=0; j<rnz; j++) {
237         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
238       }
239     }
240     *r = row; *c = col;
241   }
242   PetscFunctionReturn(0);
243 }
244 
245 PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
246 {
247   Mat_SeqBAIJ    *aa=(Mat_SeqBAIJ*)A->data;
248   const PetscInt *ai,*aj,*ajj,bs2 = aa->bs2;
249   PetscInt       bs,M,nz,idx=0,rnz,i,j,k,m;
250   PetscErrorCode ierr;
251   PetscInt       *row,*col;
252 
253   PetscFunctionBegin;
254   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
255   M = A->rmap->N/bs;
256   *v = aa->a;
257   if (reuse == MAT_INITIAL_MATRIX) {
258     ai   = aa->i; aj = aa->j;
259     nz   = bs2*aa->nz;
260     *nnz = nz;
261     ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr);
262     col  = row + nz;
263 
264     for (i=0; i<M; i++) {
265       ajj = aj + ai[i];
266       rnz = ai[i+1] - ai[i];
267       for (k=0; k<rnz; k++) {
268         for (j=0; j<bs; j++) {
269           for (m=0; m<bs; m++) {
270             row[idx]   = i*bs + m + shift;
271             col[idx++] = bs*(ajj[k]) + j + shift;
272           }
273         }
274       }
275     }
276     *r = row; *c = col;
277   }
278   PetscFunctionReturn(0);
279 }
280 
281 PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
282 {
283   const PetscInt *ai, *aj,*ajj,M=A->rmap->n;
284   PetscInt       nz,rnz,i,j;
285   PetscErrorCode ierr;
286   PetscInt       *row,*col;
287   Mat_SeqSBAIJ   *aa=(Mat_SeqSBAIJ*)A->data;
288 
289   PetscFunctionBegin;
290   *v = aa->a;
291   if (reuse == MAT_INITIAL_MATRIX) {
292     nz   = aa->nz;
293     ai   = aa->i;
294     aj   = aa->j;
295     *v   = aa->a;
296     *nnz = nz;
297     ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr);
298     col  = row + nz;
299 
300     nz = 0;
301     for (i=0; i<M; i++) {
302       rnz = ai[i+1] - ai[i];
303       ajj = aj + ai[i];
304       for (j=0; j<rnz; j++) {
305         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
306       }
307     }
308     *r = row; *c = col;
309   }
310   PetscFunctionReturn(0);
311 }
312 
313 PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
314 {
315   const PetscInt    *ai,*aj,*ajj,*adiag,M=A->rmap->n;
316   PetscInt          nz,rnz,i,j;
317   const PetscScalar *av,*v1;
318   PetscScalar       *val;
319   PetscErrorCode    ierr;
320   PetscInt          *row,*col;
321   Mat_SeqAIJ        *aa=(Mat_SeqAIJ*)A->data;
322   PetscBool         missing;
323 
324   PetscFunctionBegin;
325   ai   =aa->i; aj=aa->j;av=aa->a;
326   adiag=aa->diag;
327   ierr  = MatMissingDiagonal_SeqAIJ(A,&missing,&i);CHKERRQ(ierr);
328   if (reuse == MAT_INITIAL_MATRIX) {
329     /* count nz in the uppper triangular part of A */
330     nz = 0;
331     if (missing) {
332       for (i=0; i<M; i++) {
333         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
334           for (j=ai[i];j<ai[i+1];j++) {
335             if (aj[j] < i) continue;
336             nz++;
337           }
338         } else {
339           nz += ai[i+1] - adiag[i];
340         }
341       }
342     } else {
343       for (i=0; i<M; i++) nz += ai[i+1] - adiag[i];
344     }
345     *nnz = nz;
346 
347     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
348     col  = row + nz;
349     val  = (PetscScalar*)(col + nz);
350 
351     nz = 0;
352     if (missing) {
353       for (i=0; i<M; i++) {
354         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
355           for (j=ai[i];j<ai[i+1];j++) {
356             if (aj[j] < i) continue;
357             row[nz] = i+shift;
358             col[nz] = aj[j]+shift;
359             val[nz] = av[j];
360             nz++;
361           }
362         } else {
363           rnz = ai[i+1] - adiag[i];
364           ajj = aj + adiag[i];
365           v1  = av + adiag[i];
366           for (j=0; j<rnz; j++) {
367             row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
368           }
369         }
370       }
371     } else {
372       for (i=0; i<M; i++) {
373         rnz = ai[i+1] - adiag[i];
374         ajj = aj + adiag[i];
375         v1  = av + adiag[i];
376         for (j=0; j<rnz; j++) {
377           row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
378         }
379       }
380     }
381     *r = row; *c = col; *v = val;
382   } else {
383     nz = 0; val = *v;
384     if (missing) {
385       for (i=0; i <M; i++) {
386         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
387           for (j=ai[i];j<ai[i+1];j++) {
388             if (aj[j] < i) continue;
389             val[nz++] = av[j];
390           }
391         } else {
392           rnz = ai[i+1] - adiag[i];
393           v1  = av + adiag[i];
394           for (j=0; j<rnz; j++) {
395             val[nz++] = v1[j];
396           }
397         }
398       }
399     } else {
400       for (i=0; i <M; i++) {
401         rnz = ai[i+1] - adiag[i];
402         v1  = av + adiag[i];
403         for (j=0; j<rnz; j++) {
404           val[nz++] = v1[j];
405         }
406       }
407     }
408   }
409   PetscFunctionReturn(0);
410 }
411 
412 PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
413 {
414   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
415   PetscErrorCode    ierr;
416   PetscInt          rstart,nz,i,j,jj,irow,countA,countB;
417   PetscInt          *row,*col;
418   const PetscScalar *av, *bv,*v1,*v2;
419   PetscScalar       *val;
420   Mat_MPISBAIJ      *mat = (Mat_MPISBAIJ*)A->data;
421   Mat_SeqSBAIJ      *aa  = (Mat_SeqSBAIJ*)(mat->A)->data;
422   Mat_SeqBAIJ       *bb  = (Mat_SeqBAIJ*)(mat->B)->data;
423 
424   PetscFunctionBegin;
425   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
426   av=aa->a; bv=bb->a;
427 
428   garray = mat->garray;
429 
430   if (reuse == MAT_INITIAL_MATRIX) {
431     nz   = aa->nz + bb->nz;
432     *nnz = nz;
433     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
434     col  = row + nz;
435     val  = (PetscScalar*)(col + nz);
436 
437     *r = row; *c = col; *v = val;
438   } else {
439     row = *r; col = *c; val = *v;
440   }
441 
442   jj = 0; irow = rstart;
443   for (i=0; i<m; i++) {
444     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
445     countA = ai[i+1] - ai[i];
446     countB = bi[i+1] - bi[i];
447     bjj    = bj + bi[i];
448     v1     = av + ai[i];
449     v2     = bv + bi[i];
450 
451     /* A-part */
452     for (j=0; j<countA; j++) {
453       if (reuse == MAT_INITIAL_MATRIX) {
454         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
455       }
456       val[jj++] = v1[j];
457     }
458 
459     /* B-part */
460     for (j=0; j < countB; j++) {
461       if (reuse == MAT_INITIAL_MATRIX) {
462         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
463       }
464       val[jj++] = v2[j];
465     }
466     irow++;
467   }
468   PetscFunctionReturn(0);
469 }
470 
471 PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
472 {
473   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
474   PetscErrorCode    ierr;
475   PetscInt          rstart,nz,i,j,jj,irow,countA,countB;
476   PetscInt          *row,*col;
477   const PetscScalar *av, *bv,*v1,*v2;
478   PetscScalar       *val;
479   Mat_MPIAIJ        *mat = (Mat_MPIAIJ*)A->data;
480   Mat_SeqAIJ        *aa  = (Mat_SeqAIJ*)(mat->A)->data;
481   Mat_SeqAIJ        *bb  = (Mat_SeqAIJ*)(mat->B)->data;
482 
483   PetscFunctionBegin;
484   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
485   av=aa->a; bv=bb->a;
486 
487   garray = mat->garray;
488 
489   if (reuse == MAT_INITIAL_MATRIX) {
490     nz   = aa->nz + bb->nz;
491     *nnz = nz;
492     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
493     col  = row + nz;
494     val  = (PetscScalar*)(col + nz);
495 
496     *r = row; *c = col; *v = val;
497   } else {
498     row = *r; col = *c; val = *v;
499   }
500 
501   jj = 0; irow = rstart;
502   for (i=0; i<m; i++) {
503     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
504     countA = ai[i+1] - ai[i];
505     countB = bi[i+1] - bi[i];
506     bjj    = bj + bi[i];
507     v1     = av + ai[i];
508     v2     = bv + bi[i];
509 
510     /* A-part */
511     for (j=0; j<countA; j++) {
512       if (reuse == MAT_INITIAL_MATRIX) {
513         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
514       }
515       val[jj++] = v1[j];
516     }
517 
518     /* B-part */
519     for (j=0; j < countB; j++) {
520       if (reuse == MAT_INITIAL_MATRIX) {
521         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
522       }
523       val[jj++] = v2[j];
524     }
525     irow++;
526   }
527   PetscFunctionReturn(0);
528 }
529 
530 PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
531 {
532   Mat_MPIBAIJ       *mat    = (Mat_MPIBAIJ*)A->data;
533   Mat_SeqBAIJ       *aa     = (Mat_SeqBAIJ*)(mat->A)->data;
534   Mat_SeqBAIJ       *bb     = (Mat_SeqBAIJ*)(mat->B)->data;
535   const PetscInt    *ai     = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj;
536   const PetscInt    *garray = mat->garray,mbs=mat->mbs,rstart=A->rmap->rstart;
537   const PetscInt    bs2=mat->bs2;
538   PetscErrorCode    ierr;
539   PetscInt          bs,nz,i,j,k,n,jj,irow,countA,countB,idx;
540   PetscInt          *row,*col;
541   const PetscScalar *av=aa->a, *bv=bb->a,*v1,*v2;
542   PetscScalar       *val;
543 
544   PetscFunctionBegin;
545   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
546   if (reuse == MAT_INITIAL_MATRIX) {
547     nz   = bs2*(aa->nz + bb->nz);
548     *nnz = nz;
549     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
550     col  = row + nz;
551     val  = (PetscScalar*)(col + nz);
552 
553     *r = row; *c = col; *v = val;
554   } else {
555     row = *r; col = *c; val = *v;
556   }
557 
558   jj = 0; irow = rstart;
559   for (i=0; i<mbs; i++) {
560     countA = ai[i+1] - ai[i];
561     countB = bi[i+1] - bi[i];
562     ajj    = aj + ai[i];
563     bjj    = bj + bi[i];
564     v1     = av + bs2*ai[i];
565     v2     = bv + bs2*bi[i];
566 
567     idx = 0;
568     /* A-part */
569     for (k=0; k<countA; k++) {
570       for (j=0; j<bs; j++) {
571         for (n=0; n<bs; n++) {
572           if (reuse == MAT_INITIAL_MATRIX) {
573             row[jj] = irow + n + shift;
574             col[jj] = rstart + bs*ajj[k] + j + shift;
575           }
576           val[jj++] = v1[idx++];
577         }
578       }
579     }
580 
581     idx = 0;
582     /* B-part */
583     for (k=0; k<countB; k++) {
584       for (j=0; j<bs; j++) {
585         for (n=0; n<bs; n++) {
586           if (reuse == MAT_INITIAL_MATRIX) {
587             row[jj] = irow + n + shift;
588             col[jj] = bs*garray[bjj[k]] + j + shift;
589           }
590           val[jj++] = v2[idx++];
591         }
592       }
593     }
594     irow += bs;
595   }
596   PetscFunctionReturn(0);
597 }
598 
599 PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
600 {
601   const PetscInt    *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
602   PetscErrorCode    ierr;
603   PetscInt          rstart,nz,nza,nzb,i,j,jj,irow,countA,countB;
604   PetscInt          *row,*col;
605   const PetscScalar *av, *bv,*v1,*v2;
606   PetscScalar       *val;
607   Mat_MPIAIJ        *mat =  (Mat_MPIAIJ*)A->data;
608   Mat_SeqAIJ        *aa  =(Mat_SeqAIJ*)(mat->A)->data;
609   Mat_SeqAIJ        *bb  =(Mat_SeqAIJ*)(mat->B)->data;
610 
611   PetscFunctionBegin;
612   ai=aa->i; aj=aa->j; adiag=aa->diag;
613   bi=bb->i; bj=bb->j; garray = mat->garray;
614   av=aa->a; bv=bb->a;
615 
616   rstart = A->rmap->rstart;
617 
618   if (reuse == MAT_INITIAL_MATRIX) {
619     nza = 0;    /* num of upper triangular entries in mat->A, including diagonals */
620     nzb = 0;    /* num of upper triangular entries in mat->B */
621     for (i=0; i<m; i++) {
622       nza   += (ai[i+1] - adiag[i]);
623       countB = bi[i+1] - bi[i];
624       bjj    = bj + bi[i];
625       for (j=0; j<countB; j++) {
626         if (garray[bjj[j]] > rstart) nzb++;
627       }
628     }
629 
630     nz   = nza + nzb; /* total nz of upper triangular part of mat */
631     *nnz = nz;
632     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
633     col  = row + nz;
634     val  = (PetscScalar*)(col + nz);
635 
636     *r = row; *c = col; *v = val;
637   } else {
638     row = *r; col = *c; val = *v;
639   }
640 
641   jj = 0; irow = rstart;
642   for (i=0; i<m; i++) {
643     ajj    = aj + adiag[i];                 /* ptr to the beginning of the diagonal of this row */
644     v1     = av + adiag[i];
645     countA = ai[i+1] - adiag[i];
646     countB = bi[i+1] - bi[i];
647     bjj    = bj + bi[i];
648     v2     = bv + bi[i];
649 
650     /* A-part */
651     for (j=0; j<countA; j++) {
652       if (reuse == MAT_INITIAL_MATRIX) {
653         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
654       }
655       val[jj++] = v1[j];
656     }
657 
658     /* B-part */
659     for (j=0; j < countB; j++) {
660       if (garray[bjj[j]] > rstart) {
661         if (reuse == MAT_INITIAL_MATRIX) {
662           row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
663         }
664         val[jj++] = v2[j];
665       }
666     }
667     irow++;
668   }
669   PetscFunctionReturn(0);
670 }
671 
672 PetscErrorCode MatDestroy_MUMPS(Mat A)
673 {
674   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->data;
675   PetscErrorCode ierr;
676 
677   PetscFunctionBegin;
678   ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr);
679   ierr = VecScatterDestroy(&mumps->scat_rhs);CHKERRQ(ierr);
680   ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
681   ierr = VecDestroy(&mumps->b_seq);CHKERRQ(ierr);
682   ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr);
683   ierr = PetscFree(mumps->id.perm_in);CHKERRQ(ierr);
684   ierr = PetscFree(mumps->irn);CHKERRQ(ierr);
685   ierr = PetscFree(mumps->info);CHKERRQ(ierr);
686   ierr = MatMumpsResetSchur_Private(mumps);CHKERRQ(ierr);
687   mumps->id.job = JOB_END;
688   PetscMUMPS_c(&mumps->id);
689   ierr = MPI_Comm_free(&mumps->comm_mumps);CHKERRQ(ierr);
690   ierr = PetscFree(A->data);CHKERRQ(ierr);
691 
692   /* clear composed functions */
693   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverType_C",NULL);CHKERRQ(ierr);
694   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorSetSchurIS_C",NULL);CHKERRQ(ierr);
695   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorCreateSchurComplement_C",NULL);CHKERRQ(ierr);
696   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetIcntl_C",NULL);CHKERRQ(ierr);
697   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetIcntl_C",NULL);CHKERRQ(ierr);
698   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetCntl_C",NULL);CHKERRQ(ierr);
699   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetCntl_C",NULL);CHKERRQ(ierr);
700   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfo_C",NULL);CHKERRQ(ierr);
701   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfog_C",NULL);CHKERRQ(ierr);
702   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfo_C",NULL);CHKERRQ(ierr);
703   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfog_C",NULL);CHKERRQ(ierr);
704   PetscFunctionReturn(0);
705 }
706 
707 PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x)
708 {
709   Mat_MUMPS        *mumps=(Mat_MUMPS*)A->data;
710   PetscScalar      *array;
711   Vec              b_seq;
712   IS               is_iden,is_petsc;
713   PetscErrorCode   ierr;
714   PetscInt         i;
715   PetscBool        second_solve = PETSC_FALSE;
716   static PetscBool cite1 = PETSC_FALSE,cite2 = PETSC_FALSE;
717 
718   PetscFunctionBegin;
719   ierr = PetscCitationsRegister("@article{MUMPS01,\n  author = {P.~R. Amestoy and I.~S. Duff and J.-Y. L'Excellent and J. Koster},\n  title = {A fully asynchronous multifrontal solver using distributed dynamic scheduling},\n  journal = {SIAM Journal on Matrix Analysis and Applications},\n  volume = {23},\n  number = {1},\n  pages = {15--41},\n  year = {2001}\n}\n",&cite1);CHKERRQ(ierr);
720   ierr = PetscCitationsRegister("@article{MUMPS02,\n  author = {P.~R. Amestoy and A. Guermouche and J.-Y. L'Excellent and S. Pralet},\n  title = {Hybrid scheduling for the parallel solution of linear systems},\n  journal = {Parallel Computing},\n  volume = {32},\n  number = {2},\n  pages = {136--156},\n  year = {2006}\n}\n",&cite2);CHKERRQ(ierr);
721 
722   if (A->factorerrortype) {
723     ierr = PetscInfo2(A,"MatSolve is called with singular matrix factor, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
724     ierr = VecSetInf(x);CHKERRQ(ierr);
725     PetscFunctionReturn(0);
726   }
727 
728   mumps->id.nrhs = 1;
729   b_seq          = mumps->b_seq;
730   if (mumps->size > 1) {
731     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
732     ierr = VecScatterBegin(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
733     ierr = VecScatterEnd(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
734     if (!mumps->myid) {ierr = VecGetArray(b_seq,&array);CHKERRQ(ierr);}
735   } else {  /* size == 1 */
736     ierr = VecCopy(b,x);CHKERRQ(ierr);
737     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
738   }
739   if (!mumps->myid) { /* define rhs on the host */
740     mumps->id.nrhs = 1;
741     mumps->id.rhs = (MumpsScalar*)array;
742   }
743 
744   /*
745      handle condensation step of Schur complement (if any)
746      We set by default ICNTL(26) == -1 when Schur indices have been provided by the user.
747      According to MUMPS (5.0.0) manual, any value should be harmful during the factorization phase
748      Unless the user provides a valid value for ICNTL(26), MatSolve and MatMatSolve routines solve the full system.
749      This requires an extra call to PetscMUMPS_c and the computation of the factors for S
750   */
751   if (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2) {
752     if (mumps->size > 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Parallel Schur complements not yet supported from PETSc\n");
753     second_solve = PETSC_TRUE;
754     ierr = MatMumpsHandleSchur_Private(A,PETSC_FALSE);CHKERRQ(ierr);
755   }
756   /* solve phase */
757   /*-------------*/
758   mumps->id.job = JOB_SOLVE;
759   PetscMUMPS_c(&mumps->id);
760   if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
761 
762   /* handle expansion step of Schur complement (if any) */
763   if (second_solve) {
764     ierr = MatMumpsHandleSchur_Private(A,PETSC_TRUE);CHKERRQ(ierr);
765   }
766 
767   if (mumps->size > 1) { /* convert mumps distributed solution to petsc mpi x */
768     if (mumps->scat_sol && mumps->ICNTL9_pre != mumps->id.ICNTL(9)) {
769       /* when id.ICNTL(9) changes, the contents of lsol_loc may change (not its size, lsol_loc), recreates scat_sol */
770       ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
771     }
772     if (!mumps->scat_sol) { /* create scatter scat_sol */
773       ierr = ISCreateStride(PETSC_COMM_SELF,mumps->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */
774       for (i=0; i<mumps->id.lsol_loc; i++) {
775         mumps->id.isol_loc[i] -= 1; /* change Fortran style to C style */
776       }
777       ierr = ISCreateGeneral(PETSC_COMM_SELF,mumps->id.lsol_loc,mumps->id.isol_loc,PETSC_COPY_VALUES,&is_petsc);CHKERRQ(ierr);  /* to */
778       ierr = VecScatterCreate(mumps->x_seq,is_iden,x,is_petsc,&mumps->scat_sol);CHKERRQ(ierr);
779       ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
780       ierr = ISDestroy(&is_petsc);CHKERRQ(ierr);
781 
782       mumps->ICNTL9_pre = mumps->id.ICNTL(9); /* save current value of id.ICNTL(9) */
783     }
784 
785     ierr = VecScatterBegin(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
786     ierr = VecScatterEnd(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
787   }
788   PetscFunctionReturn(0);
789 }
790 
791 PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x)
792 {
793   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->data;
794   PetscErrorCode ierr;
795 
796   PetscFunctionBegin;
797   mumps->id.ICNTL(9) = 0;
798   ierr = MatSolve_MUMPS(A,b,x);CHKERRQ(ierr);
799   mumps->id.ICNTL(9) = 1;
800   PetscFunctionReturn(0);
801 }
802 
803 PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X)
804 {
805   PetscErrorCode ierr;
806   Mat            Bt = NULL;
807   PetscBool      flg, flgT;
808   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->data;
809   PetscInt       i,nrhs,M;
810   PetscScalar    *array,*bray;
811 
812   PetscFunctionBegin;
813   ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
814   ierr = PetscObjectTypeCompare((PetscObject)B,MATTRANSPOSEMAT,&flgT);CHKERRQ(ierr);
815   if (flgT) {
816     if (mumps->size > 1) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
817     ierr = MatTransposeGetMat(B,&Bt);CHKERRQ(ierr);
818   } else {
819     if (!flg) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
820   }
821 
822   ierr = MatGetSize(B,&M,&nrhs);CHKERRQ(ierr);
823   mumps->id.nrhs = nrhs;
824   mumps->id.lrhs = M;
825 
826   if (mumps->id.ICNTL(30)) {
827     if (mumps->size == 1 && Bt) {
828       PetscBool   done;
829       PetscScalar *aa;
830       PetscInt    spnr,*ia,*ja;
831 
832       ierr = MatSeqAIJGetArray(Bt,&aa);CHKERRQ(ierr);
833       ierr = MatGetRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&done);CHKERRQ(ierr);
834       if (!done) SETERRQ(PetscObjectComm((PetscObject)Bt),PETSC_ERR_ARG_WRONG,"Cannot get IJ structure");
835 
836       mumps->id.irhs_ptr    = ia;
837       mumps->id.irhs_sparse = ja;
838       mumps->id.nz_rhs      = ia[spnr] - 1;
839       mumps->id.rhs_sparse  = (MumpsScalar*)aa;
840       mumps->id.ICNTL(20)   = 1; /* rhs is sparse */
841 
842       /* solve phase */
843       /*-------------*/
844       mumps->id.job = JOB_SOLVE;
845       PetscMUMPS_c(&mumps->id);
846       if (mumps->id.INFOG(1) < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d INFOG(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFOG(2));
847     } else {
848       SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"not done yet");
849     }
850     PetscFunctionReturn(0);
851   } else {
852     ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
853     if (!flg) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
854   }
855 
856   if (B->rmap->n != X->rmap->n) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B and X must have same row distribution");
857 
858   if (mumps->size == 1) {
859     PetscScalar *aa;
860     PetscInt    spnr,*ia,*ja;
861     PetscBool   second_solve = PETSC_FALSE;
862 
863     /* copy B to X */
864     ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr);
865     mumps->id.rhs = (MumpsScalar*)array;
866     if (!Bt) {
867       ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr);
868       ierr = PetscMemcpy(array,bray,M*nrhs*sizeof(PetscScalar));CHKERRQ(ierr);
869       ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr);
870     } else {
871       PetscBool done;
872 
873       ierr = MatSeqAIJGetArray(Bt,&aa);CHKERRQ(ierr);
874       ierr = MatGetRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&done);CHKERRQ(ierr);
875       if (!done) SETERRQ(PetscObjectComm((PetscObject)Bt),PETSC_ERR_ARG_WRONG,"Cannot get IJ structure");
876       mumps->id.irhs_ptr    = ia;
877       mumps->id.irhs_sparse = ja;
878       mumps->id.nz_rhs      = ia[spnr] - 1;
879       mumps->id.rhs_sparse  = (MumpsScalar*)aa;
880       mumps->id.ICNTL(20)   = 1;
881     }
882     /* handle condensation step of Schur complement (if any) */
883     if (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2) {
884       second_solve = PETSC_TRUE;
885       ierr = MatMumpsHandleSchur_Private(A,PETSC_FALSE);CHKERRQ(ierr);
886     }
887     /* solve phase */
888     /*-------------*/
889     mumps->id.job = JOB_SOLVE;
890     PetscMUMPS_c(&mumps->id);
891     if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
892 
893     /* handle expansion step of Schur complement (if any) */
894     if (second_solve) {
895       ierr = MatMumpsHandleSchur_Private(A,PETSC_TRUE);CHKERRQ(ierr);
896     }
897     if (Bt) {
898       PetscBool done;
899 
900       ierr = MatSeqAIJRestoreArray(Bt,&aa);CHKERRQ(ierr);
901       ierr = MatRestoreRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&done);CHKERRQ(ierr);
902       if (!done) SETERRQ(PetscObjectComm((PetscObject)Bt),PETSC_ERR_ARG_WRONG,"Cannot restore IJ structure");
903       mumps->id.ICNTL(20) = 0;
904     }
905     ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr);
906   } else {  /*--------- parallel case --------*/
907     PetscInt       lsol_loc,nlsol_loc,*isol_loc,*idx,*iidx,*idxx,*isol_loc_save;
908     MumpsScalar    *sol_loc,*sol_loc_save;
909     IS             is_to,is_from;
910     PetscInt       k,proc,j,m;
911     const PetscInt *rstart;
912     Vec            v_mpi,b_seq,x_seq;
913     VecScatter     scat_rhs,scat_sol;
914 
915     if (mumps->size > 1 && mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Parallel Schur complements not yet supported from PETSc\n");
916 
917     /* create x_seq to hold local solution */
918     isol_loc_save = mumps->id.isol_loc; /* save it for MatSovle() */
919     sol_loc_save  = mumps->id.sol_loc;
920 
921     lsol_loc  = mumps->id.INFO(23);
922     nlsol_loc = nrhs*lsol_loc;     /* length of sol_loc */
923     ierr = PetscMalloc2(nlsol_loc,&sol_loc,nlsol_loc,&isol_loc);CHKERRQ(ierr);
924     mumps->id.sol_loc = (MumpsScalar*)sol_loc;
925     mumps->id.isol_loc = isol_loc;
926 
927     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,nlsol_loc,(PetscScalar*)sol_loc,&x_seq);CHKERRQ(ierr);
928 
929     /* copy rhs matrix B into vector v_mpi */
930     ierr = MatGetLocalSize(B,&m,NULL);CHKERRQ(ierr);
931     ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr);
932     ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)B),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi);CHKERRQ(ierr);
933     ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr);
934 
935     /* scatter v_mpi to b_seq because MUMPS only supports centralized rhs */
936     /* idx: maps from k-th index of v_mpi to (i,j)-th global entry of B;
937       iidx: inverse of idx, will be used by scattering xx_seq -> X       */
938     ierr = PetscMalloc2(nrhs*M,&idx,nrhs*M,&iidx);CHKERRQ(ierr);
939     ierr = MatGetOwnershipRanges(B,&rstart);CHKERRQ(ierr);
940     k = 0;
941     for (proc=0; proc<mumps->size; proc++){
942       for (j=0; j<nrhs; j++){
943         for (i=rstart[proc]; i<rstart[proc+1]; i++){
944           iidx[j*M + i] = k;
945           idx[k++]      = j*M + i;
946         }
947       }
948     }
949 
950     if (!mumps->myid) {
951       ierr = VecCreateSeq(PETSC_COMM_SELF,nrhs*M,&b_seq);CHKERRQ(ierr);
952       ierr = ISCreateGeneral(PETSC_COMM_SELF,nrhs*M,idx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr);
953       ierr = ISCreateStride(PETSC_COMM_SELF,nrhs*M,0,1,&is_from);CHKERRQ(ierr);
954     } else {
955       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&b_seq);CHKERRQ(ierr);
956       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_to);CHKERRQ(ierr);
957       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_from);CHKERRQ(ierr);
958     }
959     ierr = VecScatterCreate(v_mpi,is_from,b_seq,is_to,&scat_rhs);CHKERRQ(ierr);
960     ierr = VecScatterBegin(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
961     ierr = ISDestroy(&is_to);CHKERRQ(ierr);
962     ierr = ISDestroy(&is_from);CHKERRQ(ierr);
963     ierr = VecScatterEnd(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
964 
965     if (!mumps->myid) { /* define rhs on the host */
966       ierr = VecGetArray(b_seq,&bray);CHKERRQ(ierr);
967       mumps->id.rhs = (MumpsScalar*)bray;
968       ierr = VecRestoreArray(b_seq,&bray);CHKERRQ(ierr);
969     }
970 
971     /* solve phase */
972     /*-------------*/
973     mumps->id.job = JOB_SOLVE;
974     PetscMUMPS_c(&mumps->id);
975     if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
976 
977     /* scatter mumps distributed solution to petsc vector v_mpi, which shares local arrays with solution matrix X */
978     ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr);
979     ierr = VecPlaceArray(v_mpi,array);CHKERRQ(ierr);
980 
981     /* create scatter scat_sol */
982     ierr = PetscMalloc1(nlsol_loc,&idxx);CHKERRQ(ierr);
983     ierr = ISCreateStride(PETSC_COMM_SELF,nlsol_loc,0,1,&is_from);CHKERRQ(ierr);
984     for (i=0; i<lsol_loc; i++) {
985       isol_loc[i] -= 1; /* change Fortran style to C style */
986       idxx[i] = iidx[isol_loc[i]];
987       for (j=1; j<nrhs; j++){
988         idxx[j*lsol_loc+i] = iidx[isol_loc[i]+j*M];
989       }
990     }
991     ierr = ISCreateGeneral(PETSC_COMM_SELF,nlsol_loc,idxx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr);
992     ierr = VecScatterCreate(x_seq,is_from,v_mpi,is_to,&scat_sol);CHKERRQ(ierr);
993     ierr = VecScatterBegin(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
994     ierr = ISDestroy(&is_from);CHKERRQ(ierr);
995     ierr = ISDestroy(&is_to);CHKERRQ(ierr);
996     ierr = VecScatterEnd(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
997     ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr);
998 
999     /* free spaces */
1000     mumps->id.sol_loc = sol_loc_save;
1001     mumps->id.isol_loc = isol_loc_save;
1002 
1003     ierr = PetscFree2(sol_loc,isol_loc);CHKERRQ(ierr);
1004     ierr = PetscFree2(idx,iidx);CHKERRQ(ierr);
1005     ierr = PetscFree(idxx);CHKERRQ(ierr);
1006     ierr = VecDestroy(&x_seq);CHKERRQ(ierr);
1007     ierr = VecDestroy(&v_mpi);CHKERRQ(ierr);
1008     ierr = VecDestroy(&b_seq);CHKERRQ(ierr);
1009     ierr = VecScatterDestroy(&scat_rhs);CHKERRQ(ierr);
1010     ierr = VecScatterDestroy(&scat_sol);CHKERRQ(ierr);
1011   }
1012   PetscFunctionReturn(0);
1013 }
1014 
1015 #if !defined(PETSC_USE_COMPLEX)
1016 /*
1017   input:
1018    F:        numeric factor
1019   output:
1020    nneg:     total number of negative pivots
1021    nzero:    total number of zero pivots
1022    npos:     (global dimension of F) - nneg - nzero
1023 */
1024 PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
1025 {
1026   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1027   PetscErrorCode ierr;
1028   PetscMPIInt    size;
1029 
1030   PetscFunctionBegin;
1031   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)F),&size);CHKERRQ(ierr);
1032   /* 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 */
1033   if (size > 1 && mumps->id.ICNTL(13) != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"ICNTL(13)=%d. -mat_mumps_icntl_13 must be set as 1 for correct global matrix inertia\n",mumps->id.INFOG(13));
1034 
1035   if (nneg) *nneg = mumps->id.INFOG(12);
1036   if (nzero || npos) {
1037     if (mumps->id.ICNTL(24) != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"-mat_mumps_icntl_24 must be set as 1 for null pivot row detection");
1038     if (nzero) *nzero = mumps->id.INFOG(28);
1039     if (npos) *npos   = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28));
1040   }
1041   PetscFunctionReturn(0);
1042 }
1043 #endif
1044 
1045 PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
1046 {
1047   Mat_MUMPS      *mumps =(Mat_MUMPS*)(F)->data;
1048   PetscErrorCode ierr;
1049   PetscBool      isMPIAIJ;
1050 
1051   PetscFunctionBegin;
1052   if (mumps->id.INFOG(1) < 0) {
1053     if (mumps->id.INFOG(1) == -6) {
1054       ierr = PetscInfo2(A,"MatFactorNumeric is called with singular matrix structure, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
1055     }
1056     ierr = PetscInfo2(A,"MatFactorNumeric is called after analysis phase fails, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
1057     PetscFunctionReturn(0);
1058   }
1059 
1060   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1061 
1062   /* numerical factorization phase */
1063   /*-------------------------------*/
1064   mumps->id.job = JOB_FACTNUMERIC;
1065   if (!mumps->id.ICNTL(18)) { /* A is centralized */
1066     if (!mumps->myid) {
1067       mumps->id.a = (MumpsScalar*)mumps->val;
1068     }
1069   } else {
1070     mumps->id.a_loc = (MumpsScalar*)mumps->val;
1071   }
1072   PetscMUMPS_c(&mumps->id);
1073   if (mumps->id.INFOG(1) < 0) {
1074     if (A->erroriffailure) {
1075       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1076     } else {
1077       if (mumps->id.INFOG(1) == -10) { /* numerically singular matrix */
1078         ierr = PetscInfo2(F,"matrix is numerically singular, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
1079         F->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1080       } else if (mumps->id.INFOG(1) == -13) {
1081         ierr = PetscInfo2(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, cannot allocate required memory %d megabytes\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
1082         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1083       } else if (mumps->id.INFOG(1) == -8 || mumps->id.INFOG(1) == -9 || (-16 < mumps->id.INFOG(1) && mumps->id.INFOG(1) < -10) ) {
1084         ierr = PetscInfo2(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d, problem with workarray \n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
1085         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1086       } else {
1087         ierr = PetscInfo2(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
1088         F->factorerrortype = MAT_FACTOR_OTHER;
1089       }
1090     }
1091   }
1092   if (!mumps->myid && mumps->id.ICNTL(16) > 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"  mumps->id.ICNTL(16):=%d\n",mumps->id.INFOG(16));
1093 
1094   F->assembled    = PETSC_TRUE;
1095   mumps->matstruc = SAME_NONZERO_PATTERN;
1096   if (F->schur) { /* reset Schur status to unfactored */
1097     if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
1098       mumps->id.ICNTL(19) = 2;
1099       ierr = MatTranspose(F->schur,MAT_INPLACE_MATRIX,&F->schur);CHKERRQ(ierr);
1100     }
1101     ierr = MatFactorRestoreSchurComplement(F,NULL,MAT_FACTOR_SCHUR_UNFACTORED);CHKERRQ(ierr);
1102   }
1103 
1104   /* just to be sure that ICNTL(19) value returned by a call from MatMumpsGetIcntl is always consistent */
1105   if (!mumps->sym && mumps->id.ICNTL(19) && mumps->id.ICNTL(19) != 1) mumps->id.ICNTL(19) = 3;
1106 
1107   if (mumps->size > 1) {
1108     PetscInt    lsol_loc;
1109     PetscScalar *sol_loc;
1110 
1111     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr);
1112 
1113     /* distributed solution; Create x_seq=sol_loc for repeated use */
1114     if (mumps->x_seq) {
1115       ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
1116       ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr);
1117       ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr);
1118     }
1119     lsol_loc = mumps->id.INFO(23); /* length of sol_loc */
1120     ierr = PetscMalloc2(lsol_loc,&sol_loc,lsol_loc,&mumps->id.isol_loc);CHKERRQ(ierr);
1121     mumps->id.lsol_loc = lsol_loc;
1122     mumps->id.sol_loc = (MumpsScalar*)sol_loc;
1123     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq);CHKERRQ(ierr);
1124   }
1125   PetscFunctionReturn(0);
1126 }
1127 
1128 /* Sets MUMPS options from the options database */
1129 PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A)
1130 {
1131   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1132   PetscErrorCode ierr;
1133   PetscInt       icntl,info[40],i,ninfo=40;
1134   PetscBool      flg;
1135 
1136   PetscFunctionBegin;
1137   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
1138   ierr = PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);CHKERRQ(ierr);
1139   if (flg) mumps->id.ICNTL(1) = icntl;
1140   ierr = PetscOptionsInt("-mat_mumps_icntl_2","ICNTL(2): output stream for diagnostic printing, statistics, and warning","None",mumps->id.ICNTL(2),&icntl,&flg);CHKERRQ(ierr);
1141   if (flg) mumps->id.ICNTL(2) = icntl;
1142   ierr = PetscOptionsInt("-mat_mumps_icntl_3","ICNTL(3): output stream for global information, collected on the host","None",mumps->id.ICNTL(3),&icntl,&flg);CHKERRQ(ierr);
1143   if (flg) mumps->id.ICNTL(3) = icntl;
1144 
1145   ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",mumps->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
1146   if (flg) mumps->id.ICNTL(4) = icntl;
1147   if (mumps->id.ICNTL(4) || PetscLogPrintInfo) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */
1148 
1149   ierr = PetscOptionsInt("-mat_mumps_icntl_6","ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7)","None",mumps->id.ICNTL(6),&icntl,&flg);CHKERRQ(ierr);
1150   if (flg) mumps->id.ICNTL(6) = icntl;
1151 
1152   ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis","None",mumps->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
1153   if (flg) {
1154     if (icntl== 1 && mumps->size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n");
1155     else mumps->id.ICNTL(7) = icntl;
1156   }
1157 
1158   ierr = PetscOptionsInt("-mat_mumps_icntl_8","ICNTL(8): scaling strategy (-2 to 8 or 77)","None",mumps->id.ICNTL(8),&mumps->id.ICNTL(8),NULL);CHKERRQ(ierr);
1159   /* ierr = PetscOptionsInt("-mat_mumps_icntl_9","ICNTL(9): computes the solution using A or A^T","None",mumps->id.ICNTL(9),&mumps->id.ICNTL(9),NULL);CHKERRQ(ierr); handled by MatSolveTranspose_MUMPS() */
1160   ierr = PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",mumps->id.ICNTL(10),&mumps->id.ICNTL(10),NULL);CHKERRQ(ierr);
1161   ierr = PetscOptionsInt("-mat_mumps_icntl_11","ICNTL(11): statistics related to an error analysis (via -ksp_view)","None",mumps->id.ICNTL(11),&mumps->id.ICNTL(11),NULL);CHKERRQ(ierr);
1162   ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): an ordering strategy for symmetric matrices (0 to 3)","None",mumps->id.ICNTL(12),&mumps->id.ICNTL(12),NULL);CHKERRQ(ierr);
1163   ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting","None",mumps->id.ICNTL(13),&mumps->id.ICNTL(13),NULL);CHKERRQ(ierr);
1164   ierr = PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage increase in the estimated working space","None",mumps->id.ICNTL(14),&mumps->id.ICNTL(14),NULL);CHKERRQ(ierr);
1165   ierr = PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): computes the Schur complement","None",mumps->id.ICNTL(19),&mumps->id.ICNTL(19),NULL);CHKERRQ(ierr);
1166   if (mumps->id.ICNTL(19) <= 0 || mumps->id.ICNTL(19) > 3) { /* reset any schur data (if any) */
1167     ierr = MatDestroy(&F->schur);CHKERRQ(ierr);
1168     ierr = MatMumpsResetSchur_Private(mumps);CHKERRQ(ierr);
1169   }
1170   /* ierr = PetscOptionsInt("-mat_mumps_icntl_20","ICNTL(20): the format (dense or sparse) of the right-hand sides","None",mumps->id.ICNTL(20),&mumps->id.ICNTL(20),NULL);CHKERRQ(ierr); -- sparse rhs is not supported in PETSc API */
1171   /* ierr = PetscOptionsInt("-mat_mumps_icntl_21","ICNTL(21): the distribution (centralized or distributed) of the solution vectors","None",mumps->id.ICNTL(21),&mumps->id.ICNTL(21),NULL);CHKERRQ(ierr); we only use distributed solution vector */
1172 
1173   ierr = PetscOptionsInt("-mat_mumps_icntl_22","ICNTL(22): in-core/out-of-core factorization and solve (0 or 1)","None",mumps->id.ICNTL(22),&mumps->id.ICNTL(22),NULL);CHKERRQ(ierr);
1174   ierr = PetscOptionsInt("-mat_mumps_icntl_23","ICNTL(23): max size of the working memory (MB) that can allocate per processor","None",mumps->id.ICNTL(23),&mumps->id.ICNTL(23),NULL);CHKERRQ(ierr);
1175   ierr = PetscOptionsInt("-mat_mumps_icntl_24","ICNTL(24): detection of null pivot rows (0 or 1)","None",mumps->id.ICNTL(24),&mumps->id.ICNTL(24),NULL);CHKERRQ(ierr);
1176   if (mumps->id.ICNTL(24)) {
1177     mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */
1178   }
1179 
1180   ierr = PetscOptionsInt("-mat_mumps_icntl_25","ICNTL(25): computes a solution of a deficient matrix and a null space basis","None",mumps->id.ICNTL(25),&mumps->id.ICNTL(25),NULL);CHKERRQ(ierr);
1181   ierr = PetscOptionsInt("-mat_mumps_icntl_26","ICNTL(26): drives the solution phase if a Schur complement matrix","None",mumps->id.ICNTL(26),&mumps->id.ICNTL(26),NULL);CHKERRQ(ierr);
1182   ierr = PetscOptionsInt("-mat_mumps_icntl_27","ICNTL(27): the blocking size for multiple right-hand sides","None",mumps->id.ICNTL(27),&mumps->id.ICNTL(27),NULL);CHKERRQ(ierr);
1183   ierr = PetscOptionsInt("-mat_mumps_icntl_28","ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering","None",mumps->id.ICNTL(28),&mumps->id.ICNTL(28),NULL);CHKERRQ(ierr);
1184   ierr = PetscOptionsInt("-mat_mumps_icntl_29","ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis","None",mumps->id.ICNTL(29),&mumps->id.ICNTL(29),NULL);CHKERRQ(ierr);
1185   ierr = PetscOptionsInt("-mat_mumps_icntl_30","ICNTL(30): compute user-specified set of entries in inv(A)","None",mumps->id.ICNTL(30),&mumps->id.ICNTL(30),NULL);CHKERRQ(ierr);
1186   ierr = PetscOptionsInt("-mat_mumps_icntl_31","ICNTL(31): indicates which factors may be discarded during factorization","None",mumps->id.ICNTL(31),&mumps->id.ICNTL(31),NULL);CHKERRQ(ierr);
1187   /* ierr = PetscOptionsInt("-mat_mumps_icntl_32","ICNTL(32): performs the forward elemination of the right-hand sides during factorization","None",mumps->id.ICNTL(32),&mumps->id.ICNTL(32),NULL);CHKERRQ(ierr);  -- not supported by PETSc API */
1188   ierr = PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),NULL);CHKERRQ(ierr);
1189   ierr = PetscOptionsInt("-mat_mumps_icntl_35","ICNTL(35): activates Block Lock Rank (BLR) based factorization","None",mumps->id.ICNTL(35),&mumps->id.ICNTL(35),NULL);CHKERRQ(ierr);
1190 
1191   ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),NULL);CHKERRQ(ierr);
1192   ierr = PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",mumps->id.CNTL(2),&mumps->id.CNTL(2),NULL);CHKERRQ(ierr);
1193   ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),NULL);CHKERRQ(ierr);
1194   ierr = PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",mumps->id.CNTL(4),&mumps->id.CNTL(4),NULL);CHKERRQ(ierr);
1195   ierr = PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",mumps->id.CNTL(5),&mumps->id.CNTL(5),NULL);CHKERRQ(ierr);
1196   ierr = PetscOptionsReal("-mat_mumps_cntl_7","CNTL(7): dropping parameter used during BLR","None",mumps->id.CNTL(7),&mumps->id.CNTL(7),NULL);CHKERRQ(ierr);
1197 
1198   ierr = PetscOptionsString("-mat_mumps_ooc_tmpdir", "out of core directory", "None", mumps->id.ooc_tmpdir, mumps->id.ooc_tmpdir, 256, NULL);CHKERRQ(ierr);
1199 
1200   ierr = PetscOptionsIntArray("-mat_mumps_view_info","request INFO local to each processor","",info,&ninfo,NULL);CHKERRQ(ierr);
1201   if (ninfo) {
1202     if (ninfo > 40) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"number of INFO %d must <= 40\n",ninfo);
1203     ierr = PetscMalloc1(ninfo,&mumps->info);CHKERRQ(ierr);
1204     mumps->ninfo = ninfo;
1205     for (i=0; i<ninfo; i++) {
1206       if (info[i] < 0 || info[i]>40) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index of INFO %d must between 1 and 40\n",ninfo);
1207       else  mumps->info[i] = info[i];
1208     }
1209   }
1210 
1211   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1212   PetscFunctionReturn(0);
1213 }
1214 
1215 PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps)
1216 {
1217   PetscErrorCode ierr;
1218 
1219   PetscFunctionBegin;
1220   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)A), &mumps->myid);CHKERRQ(ierr);
1221   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&mumps->size);CHKERRQ(ierr);
1222   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(mumps->comm_mumps));CHKERRQ(ierr);
1223 
1224   mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps);
1225 
1226   mumps->id.job = JOB_INIT;
1227   mumps->id.par = 1;  /* host participates factorizaton and solve */
1228   mumps->id.sym = mumps->sym;
1229   PetscMUMPS_c(&mumps->id);
1230 
1231   mumps->scat_rhs     = NULL;
1232   mumps->scat_sol     = NULL;
1233 
1234   /* set PETSc-MUMPS default options - override MUMPS default */
1235   mumps->id.ICNTL(3) = 0;
1236   mumps->id.ICNTL(4) = 0;
1237   if (mumps->size == 1) {
1238     mumps->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
1239   } else {
1240     mumps->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
1241     mumps->id.ICNTL(20) = 0;   /* rhs is in dense format */
1242     mumps->id.ICNTL(21) = 1;   /* distributed solution */
1243   }
1244 
1245   /* schur */
1246   mumps->id.size_schur      = 0;
1247   mumps->id.listvar_schur   = NULL;
1248   mumps->id.schur           = NULL;
1249   mumps->sizeredrhs         = 0;
1250   mumps->schur_sol          = NULL;
1251   mumps->schur_sizesol      = 0;
1252   PetscFunctionReturn(0);
1253 }
1254 
1255 PetscErrorCode MatFactorSymbolic_MUMPS_ReportIfError(Mat F,Mat A,const MatFactorInfo *info,Mat_MUMPS *mumps)
1256 {
1257   PetscErrorCode ierr;
1258 
1259   PetscFunctionBegin;
1260   if (mumps->id.INFOG(1) < 0) {
1261     if (A->erroriffailure) {
1262       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
1263     } else {
1264       if (mumps->id.INFOG(1) == -6) {
1265         ierr = PetscInfo2(F,"matrix is singular in structure, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
1266         F->factorerrortype = MAT_FACTOR_STRUCT_ZEROPIVOT;
1267       } else if (mumps->id.INFOG(1) == -5 || mumps->id.INFOG(1) == -7) {
1268         ierr = PetscInfo2(F,"problem of workspace, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
1269         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1270       } else {
1271         ierr = PetscInfo2(F,"Error reported by MUMPS in analysis phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
1272         F->factorerrortype = MAT_FACTOR_OTHER;
1273       }
1274     }
1275   }
1276   PetscFunctionReturn(0);
1277 }
1278 
1279 /* Note Petsc r(=c) permutation is used when mumps->id.ICNTL(7)==1 with centralized assembled matrix input; otherwise r and c are ignored */
1280 PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1281 {
1282   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1283   PetscErrorCode ierr;
1284   Vec            b;
1285   IS             is_iden;
1286   const PetscInt M = A->rmap->N;
1287 
1288   PetscFunctionBegin;
1289   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1290 
1291   /* Set MUMPS options from the options database */
1292   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1293 
1294   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1295 
1296   /* analysis phase */
1297   /*----------------*/
1298   mumps->id.job = JOB_FACTSYMBOLIC;
1299   mumps->id.n   = M;
1300   switch (mumps->id.ICNTL(18)) {
1301   case 0:  /* centralized assembled matrix input */
1302     if (!mumps->myid) {
1303       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1304       if (mumps->id.ICNTL(6)>1) {
1305         mumps->id.a = (MumpsScalar*)mumps->val;
1306       }
1307       if (mumps->id.ICNTL(7) == 1) { /* use user-provide matrix ordering - assuming r = c ordering */
1308         /*
1309         PetscBool      flag;
1310         ierr = ISEqual(r,c,&flag);CHKERRQ(ierr);
1311         if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"row_perm != col_perm");
1312         ierr = ISView(r,PETSC_VIEWER_STDOUT_SELF);
1313          */
1314         if (!mumps->myid) {
1315           const PetscInt *idx;
1316           PetscInt       i,*perm_in;
1317 
1318           ierr = PetscMalloc1(M,&perm_in);CHKERRQ(ierr);
1319           ierr = ISGetIndices(r,&idx);CHKERRQ(ierr);
1320 
1321           mumps->id.perm_in = perm_in;
1322           for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */
1323           ierr = ISRestoreIndices(r,&idx);CHKERRQ(ierr);
1324         }
1325       }
1326     }
1327     break;
1328   case 3:  /* distributed assembled matrix input (size>1) */
1329     mumps->id.nz_loc = mumps->nz;
1330     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1331     if (mumps->id.ICNTL(6)>1) {
1332       mumps->id.a_loc = (MumpsScalar*)mumps->val;
1333     }
1334     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1335     if (!mumps->myid) {
1336       ierr = VecCreateSeq(PETSC_COMM_SELF,A->rmap->N,&mumps->b_seq);CHKERRQ(ierr);
1337       ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->N,0,1,&is_iden);CHKERRQ(ierr);
1338     } else {
1339       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
1340       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
1341     }
1342     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1343     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
1344     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
1345     ierr = VecDestroy(&b);CHKERRQ(ierr);
1346     break;
1347   }
1348   PetscMUMPS_c(&mumps->id);
1349   ierr = MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);CHKERRQ(ierr);
1350 
1351   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1352   F->ops->solve           = MatSolve_MUMPS;
1353   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1354   F->ops->matsolve        = MatMatSolve_MUMPS;
1355   PetscFunctionReturn(0);
1356 }
1357 
1358 /* Note the Petsc r and c permutations are ignored */
1359 PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1360 {
1361   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1362   PetscErrorCode ierr;
1363   Vec            b;
1364   IS             is_iden;
1365   const PetscInt M = A->rmap->N;
1366 
1367   PetscFunctionBegin;
1368   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1369 
1370   /* Set MUMPS options from the options database */
1371   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1372 
1373   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1374 
1375   /* analysis phase */
1376   /*----------------*/
1377   mumps->id.job = JOB_FACTSYMBOLIC;
1378   mumps->id.n   = M;
1379   switch (mumps->id.ICNTL(18)) {
1380   case 0:  /* centralized assembled matrix input */
1381     if (!mumps->myid) {
1382       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1383       if (mumps->id.ICNTL(6)>1) {
1384         mumps->id.a = (MumpsScalar*)mumps->val;
1385       }
1386     }
1387     break;
1388   case 3:  /* distributed assembled matrix input (size>1) */
1389     mumps->id.nz_loc = mumps->nz;
1390     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1391     if (mumps->id.ICNTL(6)>1) {
1392       mumps->id.a_loc = (MumpsScalar*)mumps->val;
1393     }
1394     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1395     if (!mumps->myid) {
1396       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
1397       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
1398     } else {
1399       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
1400       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
1401     }
1402     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1403     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
1404     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
1405     ierr = VecDestroy(&b);CHKERRQ(ierr);
1406     break;
1407   }
1408   PetscMUMPS_c(&mumps->id);
1409   ierr = MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);CHKERRQ(ierr);
1410 
1411   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1412   F->ops->solve           = MatSolve_MUMPS;
1413   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1414   PetscFunctionReturn(0);
1415 }
1416 
1417 /* Note the Petsc r permutation and factor info are ignored */
1418 PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
1419 {
1420   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1421   PetscErrorCode ierr;
1422   Vec            b;
1423   IS             is_iden;
1424   const PetscInt M = A->rmap->N;
1425 
1426   PetscFunctionBegin;
1427   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1428 
1429   /* Set MUMPS options from the options database */
1430   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1431 
1432   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1433 
1434   /* analysis phase */
1435   /*----------------*/
1436   mumps->id.job = JOB_FACTSYMBOLIC;
1437   mumps->id.n   = M;
1438   switch (mumps->id.ICNTL(18)) {
1439   case 0:  /* centralized assembled matrix input */
1440     if (!mumps->myid) {
1441       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1442       if (mumps->id.ICNTL(6)>1) {
1443         mumps->id.a = (MumpsScalar*)mumps->val;
1444       }
1445     }
1446     break;
1447   case 3:  /* distributed assembled matrix input (size>1) */
1448     mumps->id.nz_loc = mumps->nz;
1449     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1450     if (mumps->id.ICNTL(6)>1) {
1451       mumps->id.a_loc = (MumpsScalar*)mumps->val;
1452     }
1453     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1454     if (!mumps->myid) {
1455       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
1456       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
1457     } else {
1458       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
1459       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
1460     }
1461     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1462     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
1463     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
1464     ierr = VecDestroy(&b);CHKERRQ(ierr);
1465     break;
1466   }
1467   PetscMUMPS_c(&mumps->id);
1468   ierr = MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);CHKERRQ(ierr);
1469 
1470   F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
1471   F->ops->solve                 = MatSolve_MUMPS;
1472   F->ops->solvetranspose        = MatSolve_MUMPS;
1473   F->ops->matsolve              = MatMatSolve_MUMPS;
1474 #if defined(PETSC_USE_COMPLEX)
1475   F->ops->getinertia = NULL;
1476 #else
1477   F->ops->getinertia = MatGetInertia_SBAIJMUMPS;
1478 #endif
1479   PetscFunctionReturn(0);
1480 }
1481 
1482 PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
1483 {
1484   PetscErrorCode    ierr;
1485   PetscBool         iascii;
1486   PetscViewerFormat format;
1487   Mat_MUMPS         *mumps=(Mat_MUMPS*)A->data;
1488 
1489   PetscFunctionBegin;
1490   /* check if matrix is mumps type */
1491   if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0);
1492 
1493   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1494   if (iascii) {
1495     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1496     if (format == PETSC_VIEWER_ASCII_INFO) {
1497       ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
1498       ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                   %d \n",mumps->id.sym);CHKERRQ(ierr);
1499       ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):            %d \n",mumps->id.par);CHKERRQ(ierr);
1500       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):         %d \n",mumps->id.ICNTL(1));CHKERRQ(ierr);
1501       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg): %d \n",mumps->id.ICNTL(2));CHKERRQ(ierr);
1502       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):   %d \n",mumps->id.ICNTL(3));CHKERRQ(ierr);
1503       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):        %d \n",mumps->id.ICNTL(4));CHKERRQ(ierr);
1504       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):         %d \n",mumps->id.ICNTL(5));CHKERRQ(ierr);
1505       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):        %d \n",mumps->id.ICNTL(6));CHKERRQ(ierr);
1506       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (sequential matrix ordering):%d \n",mumps->id.ICNTL(7));CHKERRQ(ierr);
1507       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scaling strategy):        %d \n",mumps->id.ICNTL(8));CHKERRQ(ierr);
1508       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements):  %d \n",mumps->id.ICNTL(10));CHKERRQ(ierr);
1509       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):          %d \n",mumps->id.ICNTL(11));CHKERRQ(ierr);
1510       if (mumps->id.ICNTL(11)>0) {
1511         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(4) (inf norm of input mat):        %g\n",mumps->id.RINFOG(4));CHKERRQ(ierr);
1512         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(5) (inf norm of solution):         %g\n",mumps->id.RINFOG(5));CHKERRQ(ierr);
1513         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(6) (inf norm of residual):         %g\n",mumps->id.RINFOG(6));CHKERRQ(ierr);
1514         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8));CHKERRQ(ierr);
1515         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(9) (error estimate):               %g \n",mumps->id.RINFOG(9));CHKERRQ(ierr);
1516         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11));CHKERRQ(ierr);
1517       }
1518       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",mumps->id.ICNTL(12));CHKERRQ(ierr);
1519       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",mumps->id.ICNTL(13));CHKERRQ(ierr);
1520       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",mumps->id.ICNTL(14));CHKERRQ(ierr);
1521       /* ICNTL(15-17) not used */
1522       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",mumps->id.ICNTL(18));CHKERRQ(ierr);
1523       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Schur complement info):                       %d \n",mumps->id.ICNTL(19));CHKERRQ(ierr);
1524       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",mumps->id.ICNTL(20));CHKERRQ(ierr);
1525       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d \n",mumps->id.ICNTL(21));CHKERRQ(ierr);
1526       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",mumps->id.ICNTL(22));CHKERRQ(ierr);
1527       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",mumps->id.ICNTL(23));CHKERRQ(ierr);
1528 
1529       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",mumps->id.ICNTL(24));CHKERRQ(ierr);
1530       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",mumps->id.ICNTL(25));CHKERRQ(ierr);
1531       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",mumps->id.ICNTL(26));CHKERRQ(ierr);
1532       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",mumps->id.ICNTL(27));CHKERRQ(ierr);
1533       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(28) (use parallel or sequential ordering):        %d \n",mumps->id.ICNTL(28));CHKERRQ(ierr);
1534       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(29) (parallel ordering):                          %d \n",mumps->id.ICNTL(29));CHKERRQ(ierr);
1535 
1536       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(30) (user-specified set of entries in inv(A)):    %d \n",mumps->id.ICNTL(30));CHKERRQ(ierr);
1537       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(31) (factors is discarded in the solve phase):    %d \n",mumps->id.ICNTL(31));CHKERRQ(ierr);
1538       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(33) (compute determinant):                        %d \n",mumps->id.ICNTL(33));CHKERRQ(ierr);
1539       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(35) (activate BLR based factorization):           %d \n",mumps->id.ICNTL(35));CHKERRQ(ierr);
1540 
1541       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",mumps->id.CNTL(1));CHKERRQ(ierr);
1542       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2));CHKERRQ(ierr);
1543       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",mumps->id.CNTL(3));CHKERRQ(ierr);
1544       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",mumps->id.CNTL(4));CHKERRQ(ierr);
1545       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",mumps->id.CNTL(5));CHKERRQ(ierr);
1546       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(7) (dropping parameter for BLR):       %g \n",mumps->id.CNTL(7));CHKERRQ(ierr);
1547 
1548       /* infomation local to each processor */
1549       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);
1550       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
1551       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %g \n",mumps->myid,mumps->id.RINFO(1));CHKERRQ(ierr);
1552       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1553       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);
1554       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(2));CHKERRQ(ierr);
1555       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1556       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);
1557       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(3));CHKERRQ(ierr);
1558       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1559 
1560       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr);
1561       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] %d \n",mumps->myid,mumps->id.INFO(15));CHKERRQ(ierr);
1562       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1563 
1564       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr);
1565       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(16));CHKERRQ(ierr);
1566       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1567 
1568       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);
1569       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(23));CHKERRQ(ierr);
1570       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1571 
1572       if (mumps->ninfo && mumps->ninfo <= 40){
1573         PetscInt i;
1574         for (i=0; i<mumps->ninfo; i++){
1575           ierr = PetscViewerASCIIPrintf(viewer, "  INFO(%d): \n",mumps->info[i]);CHKERRQ(ierr);
1576           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(mumps->info[i]));CHKERRQ(ierr);
1577           ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1578         }
1579       }
1580 
1581 
1582       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
1583 
1584       if (!mumps->myid) { /* information from the host */
1585         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",mumps->id.RINFOG(1));CHKERRQ(ierr);
1586         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",mumps->id.RINFOG(2));CHKERRQ(ierr);
1587         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",mumps->id.RINFOG(3));CHKERRQ(ierr);
1588         ierr = PetscViewerASCIIPrintf(viewer,"  (RINFOG(12) RINFOG(13))*2^INFOG(34) (determinant): (%g,%g)*(2^%d)\n",mumps->id.RINFOG(12),mumps->id.RINFOG(13),mumps->id.INFOG(34));CHKERRQ(ierr);
1589 
1590         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(3));CHKERRQ(ierr);
1591         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(4));CHKERRQ(ierr);
1592         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",mumps->id.INFOG(5));CHKERRQ(ierr);
1593         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",mumps->id.INFOG(6));CHKERRQ(ierr);
1594         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively use after analysis): %d \n",mumps->id.INFOG(7));CHKERRQ(ierr);
1595         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",mumps->id.INFOG(8));CHKERRQ(ierr);
1596         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",mumps->id.INFOG(9));CHKERRQ(ierr);
1597         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",mumps->id.INFOG(10));CHKERRQ(ierr);
1598         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",mumps->id.INFOG(11));CHKERRQ(ierr);
1599         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",mumps->id.INFOG(12));CHKERRQ(ierr);
1600         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",mumps->id.INFOG(13));CHKERRQ(ierr);
1601         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",mumps->id.INFOG(14));CHKERRQ(ierr);
1602         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",mumps->id.INFOG(15));CHKERRQ(ierr);
1603         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",mumps->id.INFOG(16));CHKERRQ(ierr);
1604         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(17) (estimated size of all MUMPS internal data for factorization after analysis: sum over all processors): %d \n",mumps->id.INFOG(17));CHKERRQ(ierr);
1605         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(18) (size of all MUMPS internal data allocated during factorization: value on the most memory consuming processor): %d \n",mumps->id.INFOG(18));CHKERRQ(ierr);
1606         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(19) (size of all MUMPS internal data allocated during factorization: sum over all processors): %d \n",mumps->id.INFOG(19));CHKERRQ(ierr);
1607         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",mumps->id.INFOG(20));CHKERRQ(ierr);
1608         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(21) (size in MB of memory effectively used during factorization - value on the most memory consuming processor): %d \n",mumps->id.INFOG(21));CHKERRQ(ierr);
1609         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(22) (size in MB of memory effectively used during factorization - sum over all processors): %d \n",mumps->id.INFOG(22));CHKERRQ(ierr);
1610         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",mumps->id.INFOG(23));CHKERRQ(ierr);
1611         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",mumps->id.INFOG(24));CHKERRQ(ierr);
1612         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",mumps->id.INFOG(25));CHKERRQ(ierr);
1613         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(28) (after factorization: number of null pivots encountered): %d\n",mumps->id.INFOG(28));CHKERRQ(ierr);
1614         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(29) (after factorization: effective number of entries in the factors (sum over all processors)): %d\n",mumps->id.INFOG(29));CHKERRQ(ierr);
1615         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(30, 31) (after solution: size in Mbytes of memory used during solution phase): %d, %d\n",mumps->id.INFOG(30),mumps->id.INFOG(31));CHKERRQ(ierr);
1616         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(32) (after analysis: type of analysis done): %d\n",mumps->id.INFOG(32));CHKERRQ(ierr);
1617         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(33) (value used for ICNTL(8)): %d\n",mumps->id.INFOG(33));CHKERRQ(ierr);
1618         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(34) (exponent of the determinant if determinant is requested): %d\n",mumps->id.INFOG(34));CHKERRQ(ierr);
1619       }
1620     }
1621   }
1622   PetscFunctionReturn(0);
1623 }
1624 
1625 PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
1626 {
1627   Mat_MUMPS *mumps =(Mat_MUMPS*)A->data;
1628 
1629   PetscFunctionBegin;
1630   info->block_size        = 1.0;
1631   info->nz_allocated      = mumps->id.INFOG(20);
1632   info->nz_used           = mumps->id.INFOG(20);
1633   info->nz_unneeded       = 0.0;
1634   info->assemblies        = 0.0;
1635   info->mallocs           = 0.0;
1636   info->memory            = 0.0;
1637   info->fill_ratio_given  = 0;
1638   info->fill_ratio_needed = 0;
1639   info->factor_mallocs    = 0;
1640   PetscFunctionReturn(0);
1641 }
1642 
1643 /* -------------------------------------------------------------------------------------------*/
1644 PetscErrorCode MatFactorSetSchurIS_MUMPS(Mat F, IS is)
1645 {
1646   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1647   const PetscInt *idxs;
1648   PetscInt       size,i;
1649   PetscErrorCode ierr;
1650 
1651   PetscFunctionBegin;
1652   ierr = ISGetLocalSize(is,&size);CHKERRQ(ierr);
1653   if (mumps->size > 1) {
1654     PetscBool ls,gs;
1655 
1656     ls   = mumps->myid ? (size ? PETSC_FALSE : PETSC_TRUE) : PETSC_TRUE;
1657     ierr = MPI_Allreduce(&ls,&gs,1,MPIU_BOOL,MPI_LAND,mumps->comm_mumps);CHKERRQ(ierr);
1658     if (!gs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MUMPS distributed parallel Schur complements not yet supported from PETSc\n");
1659   }
1660   if (mumps->id.size_schur != size) {
1661     ierr = PetscFree2(mumps->id.listvar_schur,mumps->id.schur);CHKERRQ(ierr);
1662     mumps->id.size_schur = size;
1663     mumps->id.schur_lld  = size;
1664     ierr = PetscMalloc2(size,&mumps->id.listvar_schur,size*size,&mumps->id.schur);CHKERRQ(ierr);
1665   }
1666 
1667   /* Schur complement matrix */
1668   ierr = MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.size_schur,(PetscScalar*)mumps->id.schur,&F->schur);CHKERRQ(ierr);
1669   if (mumps->sym == 1) {
1670     ierr = MatSetOption(F->schur,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr);
1671   }
1672 
1673   /* MUMPS expects Fortran style indices */
1674   ierr = ISGetIndices(is,&idxs);CHKERRQ(ierr);
1675   ierr = PetscMemcpy(mumps->id.listvar_schur,idxs,size*sizeof(PetscInt));CHKERRQ(ierr);
1676   for (i=0;i<size;i++) mumps->id.listvar_schur[i]++;
1677   ierr = ISRestoreIndices(is,&idxs);CHKERRQ(ierr);
1678   if (mumps->size > 1) {
1679     mumps->id.ICNTL(19) = 1; /* MUMPS returns Schur centralized on the host */
1680   } else {
1681     if (F->factortype == MAT_FACTOR_LU) {
1682       mumps->id.ICNTL(19) = 3; /* MUMPS returns full matrix */
1683     } else {
1684       mumps->id.ICNTL(19) = 2; /* MUMPS returns lower triangular part */
1685     }
1686   }
1687   /* set a special value of ICNTL (not handled my MUMPS) to be used in the solve phase by PETSc */
1688   mumps->id.ICNTL(26) = -1;
1689   PetscFunctionReturn(0);
1690 }
1691 
1692 /* -------------------------------------------------------------------------------------------*/
1693 PetscErrorCode MatFactorCreateSchurComplement_MUMPS(Mat F,Mat* S)
1694 {
1695   Mat            St;
1696   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1697   PetscScalar    *array;
1698 #if defined(PETSC_USE_COMPLEX)
1699   PetscScalar    im = PetscSqrtScalar((PetscScalar)-1.0);
1700 #endif
1701   PetscErrorCode ierr;
1702 
1703   PetscFunctionBegin;
1704   if (!mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it");
1705   ierr = MatCreate(PETSC_COMM_SELF,&St);CHKERRQ(ierr);
1706   ierr = MatSetSizes(St,PETSC_DECIDE,PETSC_DECIDE,mumps->id.size_schur,mumps->id.size_schur);CHKERRQ(ierr);
1707   ierr = MatSetType(St,MATDENSE);CHKERRQ(ierr);
1708   ierr = MatSetUp(St);CHKERRQ(ierr);
1709   ierr = MatDenseGetArray(St,&array);CHKERRQ(ierr);
1710   if (!mumps->sym) { /* MUMPS always return a full matrix */
1711     if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
1712       PetscInt i,j,N=mumps->id.size_schur;
1713       for (i=0;i<N;i++) {
1714         for (j=0;j<N;j++) {
1715 #if !defined(PETSC_USE_COMPLEX)
1716           PetscScalar val = mumps->id.schur[i*N+j];
1717 #else
1718           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
1719 #endif
1720           array[j*N+i] = val;
1721         }
1722       }
1723     } else { /* stored by columns */
1724       ierr = PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
1725     }
1726   } else { /* either full or lower-triangular (not packed) */
1727     if (mumps->id.ICNTL(19) == 2) { /* lower triangular stored by columns */
1728       PetscInt i,j,N=mumps->id.size_schur;
1729       for (i=0;i<N;i++) {
1730         for (j=i;j<N;j++) {
1731 #if !defined(PETSC_USE_COMPLEX)
1732           PetscScalar val = mumps->id.schur[i*N+j];
1733 #else
1734           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
1735 #endif
1736           array[i*N+j] = val;
1737           array[j*N+i] = val;
1738         }
1739       }
1740     } else if (mumps->id.ICNTL(19) == 3) { /* full matrix */
1741       ierr = PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
1742     } else { /* ICNTL(19) == 1 lower triangular stored by rows */
1743       PetscInt i,j,N=mumps->id.size_schur;
1744       for (i=0;i<N;i++) {
1745         for (j=0;j<i+1;j++) {
1746 #if !defined(PETSC_USE_COMPLEX)
1747           PetscScalar val = mumps->id.schur[i*N+j];
1748 #else
1749           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
1750 #endif
1751           array[i*N+j] = val;
1752           array[j*N+i] = val;
1753         }
1754       }
1755     }
1756   }
1757   ierr = MatDenseRestoreArray(St,&array);CHKERRQ(ierr);
1758   *S   = St;
1759   PetscFunctionReturn(0);
1760 }
1761 
1762 /* -------------------------------------------------------------------------------------------*/
1763 PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival)
1764 {
1765   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
1766 
1767   PetscFunctionBegin;
1768   mumps->id.ICNTL(icntl) = ival;
1769   PetscFunctionReturn(0);
1770 }
1771 
1772 PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt *ival)
1773 {
1774   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
1775 
1776   PetscFunctionBegin;
1777   *ival = mumps->id.ICNTL(icntl);
1778   PetscFunctionReturn(0);
1779 }
1780 
1781 /*@
1782   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()
1783 
1784    Logically Collective on Mat
1785 
1786    Input Parameters:
1787 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1788 .  icntl - index of MUMPS parameter array ICNTL()
1789 -  ival - value of MUMPS ICNTL(icntl)
1790 
1791   Options Database:
1792 .   -mat_mumps_icntl_<icntl> <ival>
1793 
1794    Level: beginner
1795 
1796    References:
1797 .     MUMPS Users' Guide
1798 
1799 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
1800  @*/
1801 PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival)
1802 {
1803   PetscErrorCode ierr;
1804 
1805   PetscFunctionBegin;
1806   PetscValidType(F,1);
1807   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
1808   PetscValidLogicalCollectiveInt(F,icntl,2);
1809   PetscValidLogicalCollectiveInt(F,ival,3);
1810   ierr = PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr);
1811   PetscFunctionReturn(0);
1812 }
1813 
1814 /*@
1815   MatMumpsGetIcntl - Get MUMPS parameter ICNTL()
1816 
1817    Logically Collective on Mat
1818 
1819    Input Parameters:
1820 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1821 -  icntl - index of MUMPS parameter array ICNTL()
1822 
1823   Output Parameter:
1824 .  ival - value of MUMPS ICNTL(icntl)
1825 
1826    Level: beginner
1827 
1828    References:
1829 .     MUMPS Users' Guide
1830 
1831 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
1832 @*/
1833 PetscErrorCode MatMumpsGetIcntl(Mat F,PetscInt icntl,PetscInt *ival)
1834 {
1835   PetscErrorCode ierr;
1836 
1837   PetscFunctionBegin;
1838   PetscValidType(F,1);
1839   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
1840   PetscValidLogicalCollectiveInt(F,icntl,2);
1841   PetscValidIntPointer(ival,3);
1842   ierr = PetscUseMethod(F,"MatMumpsGetIcntl_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
1843   PetscFunctionReturn(0);
1844 }
1845 
1846 /* -------------------------------------------------------------------------------------------*/
1847 PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal val)
1848 {
1849   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
1850 
1851   PetscFunctionBegin;
1852   mumps->id.CNTL(icntl) = val;
1853   PetscFunctionReturn(0);
1854 }
1855 
1856 PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal *val)
1857 {
1858   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
1859 
1860   PetscFunctionBegin;
1861   *val = mumps->id.CNTL(icntl);
1862   PetscFunctionReturn(0);
1863 }
1864 
1865 /*@
1866   MatMumpsSetCntl - Set MUMPS parameter CNTL()
1867 
1868    Logically Collective on Mat
1869 
1870    Input Parameters:
1871 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1872 .  icntl - index of MUMPS parameter array CNTL()
1873 -  val - value of MUMPS CNTL(icntl)
1874 
1875   Options Database:
1876 .   -mat_mumps_cntl_<icntl> <val>
1877 
1878    Level: beginner
1879 
1880    References:
1881 .     MUMPS Users' Guide
1882 
1883 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
1884 @*/
1885 PetscErrorCode MatMumpsSetCntl(Mat F,PetscInt icntl,PetscReal val)
1886 {
1887   PetscErrorCode ierr;
1888 
1889   PetscFunctionBegin;
1890   PetscValidType(F,1);
1891   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
1892   PetscValidLogicalCollectiveInt(F,icntl,2);
1893   PetscValidLogicalCollectiveReal(F,val,3);
1894   ierr = PetscTryMethod(F,"MatMumpsSetCntl_C",(Mat,PetscInt,PetscReal),(F,icntl,val));CHKERRQ(ierr);
1895   PetscFunctionReturn(0);
1896 }
1897 
1898 /*@
1899   MatMumpsGetCntl - Get MUMPS parameter CNTL()
1900 
1901    Logically Collective on Mat
1902 
1903    Input Parameters:
1904 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1905 -  icntl - index of MUMPS parameter array CNTL()
1906 
1907   Output Parameter:
1908 .  val - value of MUMPS CNTL(icntl)
1909 
1910    Level: beginner
1911 
1912    References:
1913 .      MUMPS Users' Guide
1914 
1915 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
1916 @*/
1917 PetscErrorCode MatMumpsGetCntl(Mat F,PetscInt icntl,PetscReal *val)
1918 {
1919   PetscErrorCode ierr;
1920 
1921   PetscFunctionBegin;
1922   PetscValidType(F,1);
1923   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
1924   PetscValidLogicalCollectiveInt(F,icntl,2);
1925   PetscValidRealPointer(val,3);
1926   ierr = PetscUseMethod(F,"MatMumpsGetCntl_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
1927   PetscFunctionReturn(0);
1928 }
1929 
1930 PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F,PetscInt icntl,PetscInt *info)
1931 {
1932   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
1933 
1934   PetscFunctionBegin;
1935   *info = mumps->id.INFO(icntl);
1936   PetscFunctionReturn(0);
1937 }
1938 
1939 PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F,PetscInt icntl,PetscInt *infog)
1940 {
1941   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
1942 
1943   PetscFunctionBegin;
1944   *infog = mumps->id.INFOG(icntl);
1945   PetscFunctionReturn(0);
1946 }
1947 
1948 PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfo)
1949 {
1950   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
1951 
1952   PetscFunctionBegin;
1953   *rinfo = mumps->id.RINFO(icntl);
1954   PetscFunctionReturn(0);
1955 }
1956 
1957 PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfog)
1958 {
1959   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
1960 
1961   PetscFunctionBegin;
1962   *rinfog = mumps->id.RINFOG(icntl);
1963   PetscFunctionReturn(0);
1964 }
1965 
1966 /*@
1967   MatMumpsGetInfo - Get MUMPS parameter INFO()
1968 
1969    Logically Collective on Mat
1970 
1971    Input Parameters:
1972 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1973 -  icntl - index of MUMPS parameter array INFO()
1974 
1975   Output Parameter:
1976 .  ival - value of MUMPS INFO(icntl)
1977 
1978    Level: beginner
1979 
1980    References:
1981 .      MUMPS Users' Guide
1982 
1983 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
1984 @*/
1985 PetscErrorCode MatMumpsGetInfo(Mat F,PetscInt icntl,PetscInt *ival)
1986 {
1987   PetscErrorCode ierr;
1988 
1989   PetscFunctionBegin;
1990   PetscValidType(F,1);
1991   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
1992   PetscValidIntPointer(ival,3);
1993   ierr = PetscUseMethod(F,"MatMumpsGetInfo_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
1994   PetscFunctionReturn(0);
1995 }
1996 
1997 /*@
1998   MatMumpsGetInfog - Get MUMPS parameter INFOG()
1999 
2000    Logically Collective on Mat
2001 
2002    Input Parameters:
2003 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2004 -  icntl - index of MUMPS parameter array INFOG()
2005 
2006   Output Parameter:
2007 .  ival - value of MUMPS INFOG(icntl)
2008 
2009    Level: beginner
2010 
2011    References:
2012 .      MUMPS Users' Guide
2013 
2014 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2015 @*/
2016 PetscErrorCode MatMumpsGetInfog(Mat F,PetscInt icntl,PetscInt *ival)
2017 {
2018   PetscErrorCode ierr;
2019 
2020   PetscFunctionBegin;
2021   PetscValidType(F,1);
2022   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2023   PetscValidIntPointer(ival,3);
2024   ierr = PetscUseMethod(F,"MatMumpsGetInfog_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
2025   PetscFunctionReturn(0);
2026 }
2027 
2028 /*@
2029   MatMumpsGetRinfo - Get MUMPS parameter RINFO()
2030 
2031    Logically Collective on Mat
2032 
2033    Input Parameters:
2034 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2035 -  icntl - index of MUMPS parameter array RINFO()
2036 
2037   Output Parameter:
2038 .  val - value of MUMPS RINFO(icntl)
2039 
2040    Level: beginner
2041 
2042    References:
2043 .       MUMPS Users' Guide
2044 
2045 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2046 @*/
2047 PetscErrorCode MatMumpsGetRinfo(Mat F,PetscInt icntl,PetscReal *val)
2048 {
2049   PetscErrorCode ierr;
2050 
2051   PetscFunctionBegin;
2052   PetscValidType(F,1);
2053   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2054   PetscValidRealPointer(val,3);
2055   ierr = PetscUseMethod(F,"MatMumpsGetRinfo_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
2056   PetscFunctionReturn(0);
2057 }
2058 
2059 /*@
2060   MatMumpsGetRinfog - Get MUMPS parameter RINFOG()
2061 
2062    Logically Collective on Mat
2063 
2064    Input Parameters:
2065 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2066 -  icntl - index of MUMPS parameter array RINFOG()
2067 
2068   Output Parameter:
2069 .  val - value of MUMPS RINFOG(icntl)
2070 
2071    Level: beginner
2072 
2073    References:
2074 .      MUMPS Users' Guide
2075 
2076 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2077 @*/
2078 PetscErrorCode MatMumpsGetRinfog(Mat F,PetscInt icntl,PetscReal *val)
2079 {
2080   PetscErrorCode ierr;
2081 
2082   PetscFunctionBegin;
2083   PetscValidType(F,1);
2084   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2085   PetscValidRealPointer(val,3);
2086   ierr = PetscUseMethod(F,"MatMumpsGetRinfog_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
2087   PetscFunctionReturn(0);
2088 }
2089 
2090 /*MC
2091   MATSOLVERMUMPS -  A matrix type providing direct solvers (LU and Cholesky) for
2092   distributed and sequential matrices via the external package MUMPS.
2093 
2094   Works with MATAIJ and MATSBAIJ matrices
2095 
2096   Use ./configure --download-mumps --download-scalapack --download-parmetis --download-metis --download-ptscotch  to have PETSc installed with MUMPS
2097 
2098   Use -pc_type cholesky or lu -pc_factor_mat_solver_type mumps to use this direct solver
2099 
2100   Options Database Keys:
2101 +  -mat_mumps_icntl_1 - ICNTL(1): output stream for error messages
2102 .  -mat_mumps_icntl_2 - ICNTL(2): output stream for diagnostic printing, statistics, and warning
2103 .  -mat_mumps_icntl_3 -  ICNTL(3): output stream for global information, collected on the host
2104 .  -mat_mumps_icntl_4 -  ICNTL(4): level of printing (0 to 4)
2105 .  -mat_mumps_icntl_6 - ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7)
2106 .  -mat_mumps_icntl_7 - ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis
2107 .  -mat_mumps_icntl_8  - ICNTL(8): scaling strategy (-2 to 8 or 77)
2108 .  -mat_mumps_icntl_10  - ICNTL(10): max num of refinements
2109 .  -mat_mumps_icntl_11  - ICNTL(11): statistics related to an error analysis (via -ksp_view)
2110 .  -mat_mumps_icntl_12  - ICNTL(12): an ordering strategy for symmetric matrices (0 to 3)
2111 .  -mat_mumps_icntl_13  - ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting
2112 .  -mat_mumps_icntl_14  - ICNTL(14): percentage increase in the estimated working space
2113 .  -mat_mumps_icntl_19  - ICNTL(19): computes the Schur complement
2114 .  -mat_mumps_icntl_22  - ICNTL(22): in-core/out-of-core factorization and solve (0 or 1)
2115 .  -mat_mumps_icntl_23  - ICNTL(23): max size of the working memory (MB) that can allocate per processor
2116 .  -mat_mumps_icntl_24  - ICNTL(24): detection of null pivot rows (0 or 1)
2117 .  -mat_mumps_icntl_25  - ICNTL(25): compute a solution of a deficient matrix and a null space basis
2118 .  -mat_mumps_icntl_26  - ICNTL(26): drives the solution phase if a Schur complement matrix
2119 .  -mat_mumps_icntl_28  - ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering
2120 .  -mat_mumps_icntl_29 - ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis
2121 .  -mat_mumps_icntl_30 - ICNTL(30): compute user-specified set of entries in inv(A)
2122 .  -mat_mumps_icntl_31 - ICNTL(31): indicates which factors may be discarded during factorization
2123 .  -mat_mumps_icntl_33 - ICNTL(33): compute determinant
2124 .  -mat_mumps_cntl_1  - CNTL(1): relative pivoting threshold
2125 .  -mat_mumps_cntl_2  -  CNTL(2): stopping criterion of refinement
2126 .  -mat_mumps_cntl_3 - CNTL(3): absolute pivoting threshold
2127 .  -mat_mumps_cntl_4 - CNTL(4): value for static pivoting
2128 -  -mat_mumps_cntl_5 - CNTL(5): fixation for null pivots
2129 
2130   Level: beginner
2131 
2132     Notes: When a MUMPS factorization fails inside a KSP solve, for example with a KSP_DIVERGED_PCSETUP_FAILED, one can find the MUMPS information about the failure by calling
2133 $          KSPGetPC(ksp,&pc);
2134 $          PCFactorGetMatrix(pc,&mat);
2135 $          MatMumpsGetInfo(mat,....);
2136 $          MatMumpsGetInfog(mat,....); etc.
2137            Or you can run with -ksp_error_if_not_converged and the program will be stopped and the information printed in the error message.
2138 
2139 .seealso: PCFactorSetMatSolverType(), MatSolverType, MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog(), KSPGetPC(), PCGetFactor(), PCFactorGetMatrix()
2140 
2141 M*/
2142 
2143 static PetscErrorCode MatFactorGetSolverType_mumps(Mat A,MatSolverType *type)
2144 {
2145   PetscFunctionBegin;
2146   *type = MATSOLVERMUMPS;
2147   PetscFunctionReturn(0);
2148 }
2149 
2150 /* MatGetFactor for Seq and MPI AIJ matrices */
2151 static PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
2152 {
2153   Mat            B;
2154   PetscErrorCode ierr;
2155   Mat_MUMPS      *mumps;
2156   PetscBool      isSeqAIJ;
2157 
2158   PetscFunctionBegin;
2159   /* Create the factorization matrix */
2160   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr);
2161   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
2162   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
2163   ierr = PetscStrallocpy("mumps",&((PetscObject)B)->type_name);CHKERRQ(ierr);
2164   ierr = MatSetUp(B);CHKERRQ(ierr);
2165 
2166   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
2167 
2168   B->ops->view        = MatView_MUMPS;
2169   B->ops->getinfo     = MatGetInfo_MUMPS;
2170 
2171   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);CHKERRQ(ierr);
2172   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr);
2173   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr);
2174   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
2175   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
2176   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
2177   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
2178   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
2179   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
2180   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
2181   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
2182 
2183   if (ftype == MAT_FACTOR_LU) {
2184     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
2185     B->factortype            = MAT_FACTOR_LU;
2186     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
2187     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
2188     mumps->sym = 0;
2189   } else {
2190     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2191     B->factortype                  = MAT_FACTOR_CHOLESKY;
2192     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
2193     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
2194 #if defined(PETSC_USE_COMPLEX)
2195     mumps->sym = 2;
2196 #else
2197     if (A->spd_set && A->spd) mumps->sym = 1;
2198     else                      mumps->sym = 2;
2199 #endif
2200   }
2201 
2202   /* set solvertype */
2203   ierr = PetscFree(B->solvertype);CHKERRQ(ierr);
2204   ierr = PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);CHKERRQ(ierr);
2205 
2206   mumps->isAIJ    = PETSC_TRUE;
2207   B->ops->destroy = MatDestroy_MUMPS;
2208   B->data        = (void*)mumps;
2209 
2210   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
2211 
2212   *F = B;
2213   PetscFunctionReturn(0);
2214 }
2215 
2216 /* MatGetFactor for Seq and MPI SBAIJ matrices */
2217 static PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
2218 {
2219   Mat            B;
2220   PetscErrorCode ierr;
2221   Mat_MUMPS      *mumps;
2222   PetscBool      isSeqSBAIJ;
2223 
2224   PetscFunctionBegin;
2225   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
2226   if (A->rmap->bs > 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with block size > 1 with MUMPS Cholesky, use AIJ matrix instead");
2227   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
2228   /* Create the factorization matrix */
2229   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
2230   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
2231   ierr = PetscStrallocpy("mumps",&((PetscObject)B)->type_name);CHKERRQ(ierr);
2232   ierr = MatSetUp(B);CHKERRQ(ierr);
2233 
2234   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
2235   if (isSeqSBAIJ) {
2236     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
2237   } else {
2238     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
2239   }
2240 
2241   B->ops->getinfo                = MatGetInfo_External;
2242   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2243   B->ops->view                   = MatView_MUMPS;
2244 
2245   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);CHKERRQ(ierr);
2246   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr);
2247   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr);
2248   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
2249   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
2250   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
2251   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
2252   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
2253   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
2254   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
2255   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
2256 
2257   B->factortype = MAT_FACTOR_CHOLESKY;
2258 #if defined(PETSC_USE_COMPLEX)
2259   mumps->sym = 2;
2260 #else
2261   if (A->spd_set && A->spd) mumps->sym = 1;
2262   else                      mumps->sym = 2;
2263 #endif
2264 
2265   /* set solvertype */
2266   ierr = PetscFree(B->solvertype);CHKERRQ(ierr);
2267   ierr = PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);CHKERRQ(ierr);
2268 
2269   mumps->isAIJ    = PETSC_FALSE;
2270   B->ops->destroy = MatDestroy_MUMPS;
2271   B->data        = (void*)mumps;
2272 
2273   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
2274 
2275   *F = B;
2276   PetscFunctionReturn(0);
2277 }
2278 
2279 static PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
2280 {
2281   Mat            B;
2282   PetscErrorCode ierr;
2283   Mat_MUMPS      *mumps;
2284   PetscBool      isSeqBAIJ;
2285 
2286   PetscFunctionBegin;
2287   /* Create the factorization matrix */
2288   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr);
2289   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
2290   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
2291   ierr = PetscStrallocpy("mumps",&((PetscObject)B)->type_name);CHKERRQ(ierr);
2292   ierr = MatSetUp(B);CHKERRQ(ierr);
2293 
2294   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
2295   if (ftype == MAT_FACTOR_LU) {
2296     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
2297     B->factortype            = MAT_FACTOR_LU;
2298     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
2299     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
2300     mumps->sym = 0;
2301   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n");
2302 
2303   B->ops->getinfo     = MatGetInfo_External;
2304   B->ops->view        = MatView_MUMPS;
2305 
2306   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);CHKERRQ(ierr);
2307   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr);
2308   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr);
2309   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
2310   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
2311   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
2312   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
2313   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
2314   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
2315   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
2316   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
2317 
2318   /* set solvertype */
2319   ierr = PetscFree(B->solvertype);CHKERRQ(ierr);
2320   ierr = PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);CHKERRQ(ierr);
2321 
2322   mumps->isAIJ    = PETSC_TRUE;
2323   B->ops->destroy = MatDestroy_MUMPS;
2324   B->data        = (void*)mumps;
2325 
2326   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
2327 
2328   *F = B;
2329   PetscFunctionReturn(0);
2330 }
2331 
2332 PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_MUMPS(void)
2333 {
2334   PetscErrorCode ierr;
2335 
2336   PetscFunctionBegin;
2337   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr);
2338   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr);
2339   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr);
2340   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr);
2341   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATMPISBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr);
2342   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr);
2343   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr);
2344   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr);
2345   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr);
2346   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr);
2347   PetscFunctionReturn(0);
2348 }
2349 
2350