xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision 87b22cf42b04ab2b225e98c30ae28067ddd8708c)
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   if (mumps->id.ICNTL(30)) {
823     if (mumps->size == 1 && Bt) {
824       PetscBool   done;
825       PetscScalar *aa;
826       PetscInt    spnr,*ia,*ja;
827 
828       ierr = MatSeqAIJGetArray(Bt,&aa);CHKERRQ(ierr);
829       ierr = MatGetRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&done);CHKERRQ(ierr);
830       if (!done) SETERRQ(PetscObjectComm((PetscObject)Bt),PETSC_ERR_ARG_WRONG,"Cannot get IJ structure");
831       mumps->id.irhs_ptr    = ia;
832       mumps->id.irhs_sparse = ja;
833       mumps->id.nz_rhs      = ia[spnr] - 1;
834       mumps->id.rhs_sparse  = (MumpsScalar*)aa;
835       mumps->id.ICNTL(20)   = 1; /* rhs is sparse */
836 
837       /* solve phase */
838       /*-------------*/
839       mumps->id.job = JOB_SOLVE;
840       PetscMUMPS_c(&mumps->id);
841       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));
842 
843       printf("not done yet\n");
844     } else {
845       SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"not done yet");
846     }
847     PetscFunctionReturn(0);
848   } else {
849     ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
850     if (!flg) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
851   }
852 
853   if (B->rmap->n != X->rmap->n) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B and X must have same row distribution");
854 
855   ierr = MatGetSize(B,&M,&nrhs);CHKERRQ(ierr);
856   mumps->id.nrhs = nrhs;
857   mumps->id.lrhs = M;
858 
859   if (mumps->size == 1) {
860     PetscScalar *aa;
861     PetscInt    spnr,*ia,*ja;
862     PetscBool   second_solve = PETSC_FALSE;
863 
864     /* copy B to X */
865     ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr);
866     mumps->id.rhs = (MumpsScalar*)array;
867     if (!Bt) {
868       ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr);
869       ierr = PetscMemcpy(array,bray,M*nrhs*sizeof(PetscScalar));CHKERRQ(ierr);
870       ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr);
871     } else {
872       PetscBool done;
873 
874       ierr = MatSeqAIJGetArray(Bt,&aa);CHKERRQ(ierr);
875       ierr = MatGetRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&done);CHKERRQ(ierr);
876       if (!done) SETERRQ(PetscObjectComm((PetscObject)Bt),PETSC_ERR_ARG_WRONG,"Cannot get IJ structure");
877       mumps->id.irhs_ptr    = ia;
878       mumps->id.irhs_sparse = ja;
879       mumps->id.nz_rhs      = ia[spnr] - 1;
880       mumps->id.rhs_sparse  = (MumpsScalar*)aa;
881       mumps->id.ICNTL(20)   = 1;
882     }
883     /* handle condensation step of Schur complement (if any) */
884     if (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2) {
885       second_solve = PETSC_TRUE;
886       ierr = MatMumpsHandleSchur_Private(A,PETSC_FALSE);CHKERRQ(ierr);
887     }
888     /* solve phase */
889     /*-------------*/
890     mumps->id.job = JOB_SOLVE;
891     PetscMUMPS_c(&mumps->id);
892     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));
893 
894     /* handle expansion step of Schur complement (if any) */
895     if (second_solve) {
896       ierr = MatMumpsHandleSchur_Private(A,PETSC_TRUE);CHKERRQ(ierr);
897     }
898     if (Bt) {
899       PetscBool done;
900 
901       ierr = MatSeqAIJRestoreArray(Bt,&aa);CHKERRQ(ierr);
902       ierr = MatRestoreRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&done);CHKERRQ(ierr);
903       if (!done) SETERRQ(PetscObjectComm((PetscObject)Bt),PETSC_ERR_ARG_WRONG,"Cannot restore IJ structure");
904       mumps->id.ICNTL(20) = 0;
905     }
906     ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr);
907   } else {  /*--------- parallel case --------*/
908     PetscInt       lsol_loc,nlsol_loc,*isol_loc,*idx,*iidx,*idxx,*isol_loc_save;
909     MumpsScalar    *sol_loc,*sol_loc_save;
910     IS             is_to,is_from;
911     PetscInt       k,proc,j,m;
912     const PetscInt *rstart;
913     Vec            v_mpi,b_seq,x_seq;
914     VecScatter     scat_rhs,scat_sol;
915 
916     if (mumps->size > 1 && mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Parallel Schur complements not yet supported from PETSc\n");
917 
918     /* create x_seq to hold local solution */
919     isol_loc_save = mumps->id.isol_loc; /* save it for MatSovle() */
920     sol_loc_save  = mumps->id.sol_loc;
921 
922     lsol_loc  = mumps->id.INFO(23);
923     nlsol_loc = nrhs*lsol_loc;     /* length of sol_loc */
924     ierr = PetscMalloc2(nlsol_loc,&sol_loc,nlsol_loc,&isol_loc);CHKERRQ(ierr);
925     mumps->id.sol_loc = (MumpsScalar*)sol_loc;
926     mumps->id.isol_loc = isol_loc;
927 
928     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,nlsol_loc,(PetscScalar*)sol_loc,&x_seq);CHKERRQ(ierr);
929 
930     /* copy rhs matrix B into vector v_mpi */
931     ierr = MatGetLocalSize(B,&m,NULL);CHKERRQ(ierr);
932     ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr);
933     ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)B),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi);CHKERRQ(ierr);
934     ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr);
935 
936     /* scatter v_mpi to b_seq because MUMPS only supports centralized rhs */
937     /* idx: maps from k-th index of v_mpi to (i,j)-th global entry of B;
938       iidx: inverse of idx, will be used by scattering xx_seq -> X       */
939     ierr = PetscMalloc2(nrhs*M,&idx,nrhs*M,&iidx);CHKERRQ(ierr);
940     ierr = MatGetOwnershipRanges(B,&rstart);CHKERRQ(ierr);
941     k = 0;
942     for (proc=0; proc<mumps->size; proc++){
943       for (j=0; j<nrhs; j++){
944         for (i=rstart[proc]; i<rstart[proc+1]; i++){
945           iidx[j*M + i] = k;
946           idx[k++]      = j*M + i;
947         }
948       }
949     }
950 
951     if (!mumps->myid) {
952       ierr = VecCreateSeq(PETSC_COMM_SELF,nrhs*M,&b_seq);CHKERRQ(ierr);
953       ierr = ISCreateGeneral(PETSC_COMM_SELF,nrhs*M,idx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr);
954       ierr = ISCreateStride(PETSC_COMM_SELF,nrhs*M,0,1,&is_from);CHKERRQ(ierr);
955     } else {
956       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&b_seq);CHKERRQ(ierr);
957       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_to);CHKERRQ(ierr);
958       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_from);CHKERRQ(ierr);
959     }
960     ierr = VecScatterCreate(v_mpi,is_from,b_seq,is_to,&scat_rhs);CHKERRQ(ierr);
961     ierr = VecScatterBegin(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
962     ierr = ISDestroy(&is_to);CHKERRQ(ierr);
963     ierr = ISDestroy(&is_from);CHKERRQ(ierr);
964     ierr = VecScatterEnd(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
965 
966     if (!mumps->myid) { /* define rhs on the host */
967       ierr = VecGetArray(b_seq,&bray);CHKERRQ(ierr);
968       mumps->id.rhs = (MumpsScalar*)bray;
969       ierr = VecRestoreArray(b_seq,&bray);CHKERRQ(ierr);
970     }
971 
972     /* solve phase */
973     /*-------------*/
974     mumps->id.job = JOB_SOLVE;
975     PetscMUMPS_c(&mumps->id);
976     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));
977 
978     /* scatter mumps distributed solution to petsc vector v_mpi, which shares local arrays with solution matrix X */
979     ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr);
980     ierr = VecPlaceArray(v_mpi,array);CHKERRQ(ierr);
981 
982     /* create scatter scat_sol */
983     ierr = PetscMalloc1(nlsol_loc,&idxx);CHKERRQ(ierr);
984     ierr = ISCreateStride(PETSC_COMM_SELF,nlsol_loc,0,1,&is_from);CHKERRQ(ierr);
985     for (i=0; i<lsol_loc; i++) {
986       isol_loc[i] -= 1; /* change Fortran style to C style */
987       idxx[i] = iidx[isol_loc[i]];
988       for (j=1; j<nrhs; j++){
989         idxx[j*lsol_loc+i] = iidx[isol_loc[i]+j*M];
990       }
991     }
992     ierr = ISCreateGeneral(PETSC_COMM_SELF,nlsol_loc,idxx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr);
993     ierr = VecScatterCreate(x_seq,is_from,v_mpi,is_to,&scat_sol);CHKERRQ(ierr);
994     ierr = VecScatterBegin(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
995     ierr = ISDestroy(&is_from);CHKERRQ(ierr);
996     ierr = ISDestroy(&is_to);CHKERRQ(ierr);
997     ierr = VecScatterEnd(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
998     ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr);
999 
1000     /* free spaces */
1001     mumps->id.sol_loc = sol_loc_save;
1002     mumps->id.isol_loc = isol_loc_save;
1003 
1004     ierr = PetscFree2(sol_loc,isol_loc);CHKERRQ(ierr);
1005     ierr = PetscFree2(idx,iidx);CHKERRQ(ierr);
1006     ierr = PetscFree(idxx);CHKERRQ(ierr);
1007     ierr = VecDestroy(&x_seq);CHKERRQ(ierr);
1008     ierr = VecDestroy(&v_mpi);CHKERRQ(ierr);
1009     ierr = VecDestroy(&b_seq);CHKERRQ(ierr);
1010     ierr = VecScatterDestroy(&scat_rhs);CHKERRQ(ierr);
1011     ierr = VecScatterDestroy(&scat_sol);CHKERRQ(ierr);
1012   }
1013   PetscFunctionReturn(0);
1014 }
1015 
1016 #if !defined(PETSC_USE_COMPLEX)
1017 /*
1018   input:
1019    F:        numeric factor
1020   output:
1021    nneg:     total number of negative pivots
1022    nzero:    total number of zero pivots
1023    npos:     (global dimension of F) - nneg - nzero
1024 */
1025 PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
1026 {
1027   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1028   PetscErrorCode ierr;
1029   PetscMPIInt    size;
1030 
1031   PetscFunctionBegin;
1032   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)F),&size);CHKERRQ(ierr);
1033   /* 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 */
1034   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));
1035 
1036   if (nneg) *nneg = mumps->id.INFOG(12);
1037   if (nzero || npos) {
1038     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");
1039     if (nzero) *nzero = mumps->id.INFOG(28);
1040     if (npos) *npos   = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28));
1041   }
1042   PetscFunctionReturn(0);
1043 }
1044 #endif
1045 
1046 PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
1047 {
1048   Mat_MUMPS      *mumps =(Mat_MUMPS*)(F)->data;
1049   PetscErrorCode ierr;
1050   PetscBool      isMPIAIJ;
1051 
1052   PetscFunctionBegin;
1053   if (mumps->id.INFOG(1) < 0) {
1054     if (mumps->id.INFOG(1) == -6) {
1055       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);
1056     }
1057     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);
1058     PetscFunctionReturn(0);
1059   }
1060 
1061   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1062 
1063   /* numerical factorization phase */
1064   /*-------------------------------*/
1065   mumps->id.job = JOB_FACTNUMERIC;
1066   if (!mumps->id.ICNTL(18)) { /* A is centralized */
1067     if (!mumps->myid) {
1068       mumps->id.a = (MumpsScalar*)mumps->val;
1069     }
1070   } else {
1071     mumps->id.a_loc = (MumpsScalar*)mumps->val;
1072   }
1073   PetscMUMPS_c(&mumps->id);
1074   if (mumps->id.INFOG(1) < 0) {
1075     if (A->erroriffailure) {
1076       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));
1077     } else {
1078       if (mumps->id.INFOG(1) == -10) { /* numerically singular matrix */
1079         ierr = PetscInfo2(F,"matrix is numerically singular, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
1080         F->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1081       } else if (mumps->id.INFOG(1) == -13) {
1082         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);
1083         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1084       } else if (mumps->id.INFOG(1) == -8 || mumps->id.INFOG(1) == -9 || (-16 < mumps->id.INFOG(1) && mumps->id.INFOG(1) < -10) ) {
1085         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);
1086         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1087       } else {
1088         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);
1089         F->factorerrortype = MAT_FACTOR_OTHER;
1090       }
1091     }
1092   }
1093   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));
1094 
1095   F->assembled    = PETSC_TRUE;
1096   mumps->matstruc = SAME_NONZERO_PATTERN;
1097   if (F->schur) { /* reset Schur status to unfactored */
1098     if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
1099       mumps->id.ICNTL(19) = 2;
1100       ierr = MatTranspose(F->schur,MAT_INPLACE_MATRIX,&F->schur);CHKERRQ(ierr);
1101     }
1102     ierr = MatFactorRestoreSchurComplement(F,NULL,MAT_FACTOR_SCHUR_UNFACTORED);CHKERRQ(ierr);
1103   }
1104 
1105   /* just to be sure that ICNTL(19) value returned by a call from MatMumpsGetIcntl is always consistent */
1106   if (!mumps->sym && mumps->id.ICNTL(19) && mumps->id.ICNTL(19) != 1) mumps->id.ICNTL(19) = 3;
1107 
1108   if (mumps->size > 1) {
1109     PetscInt    lsol_loc;
1110     PetscScalar *sol_loc;
1111 
1112     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr);
1113 
1114     /* distributed solution; Create x_seq=sol_loc for repeated use */
1115     if (mumps->x_seq) {
1116       ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
1117       ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr);
1118       ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr);
1119     }
1120     lsol_loc = mumps->id.INFO(23); /* length of sol_loc */
1121     ierr = PetscMalloc2(lsol_loc,&sol_loc,lsol_loc,&mumps->id.isol_loc);CHKERRQ(ierr);
1122     mumps->id.lsol_loc = lsol_loc;
1123     mumps->id.sol_loc = (MumpsScalar*)sol_loc;
1124     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq);CHKERRQ(ierr);
1125   }
1126   PetscFunctionReturn(0);
1127 }
1128 
1129 /* Sets MUMPS options from the options database */
1130 PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A)
1131 {
1132   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1133   PetscErrorCode ierr;
1134   PetscInt       icntl,info[40],i,ninfo=40;
1135   PetscBool      flg;
1136 
1137   PetscFunctionBegin;
1138   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
1139   ierr = PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);CHKERRQ(ierr);
1140   if (flg) mumps->id.ICNTL(1) = icntl;
1141   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);
1142   if (flg) mumps->id.ICNTL(2) = icntl;
1143   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);
1144   if (flg) mumps->id.ICNTL(3) = icntl;
1145 
1146   ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",mumps->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
1147   if (flg) mumps->id.ICNTL(4) = icntl;
1148   if (mumps->id.ICNTL(4) || PetscLogPrintInfo) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */
1149 
1150   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);
1151   if (flg) mumps->id.ICNTL(6) = icntl;
1152 
1153   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);
1154   if (flg) {
1155     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");
1156     else mumps->id.ICNTL(7) = icntl;
1157   }
1158 
1159   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);
1160   /* 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() */
1161   ierr = PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",mumps->id.ICNTL(10),&mumps->id.ICNTL(10),NULL);CHKERRQ(ierr);
1162   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);
1163   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);
1164   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);
1165   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);
1166   ierr = PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): computes the Schur complement","None",mumps->id.ICNTL(19),&mumps->id.ICNTL(19),NULL);CHKERRQ(ierr);
1167   if (mumps->id.ICNTL(19) <= 0 || mumps->id.ICNTL(19) > 3) { /* reset any schur data (if any) */
1168     ierr = MatDestroy(&F->schur);CHKERRQ(ierr);
1169     ierr = MatMumpsResetSchur_Private(mumps);CHKERRQ(ierr);
1170   }
1171   /* 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 */
1172   /* 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 */
1173 
1174   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);
1175   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);
1176   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);
1177   if (mumps->id.ICNTL(24)) {
1178     mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */
1179   }
1180 
1181   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);
1182   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);
1183   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);
1184   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);
1185   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);
1186   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);
1187   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);
1188   /* 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 */
1189   ierr = PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),NULL);CHKERRQ(ierr);
1190   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);
1191 
1192   ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),NULL);CHKERRQ(ierr);
1193   ierr = PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",mumps->id.CNTL(2),&mumps->id.CNTL(2),NULL);CHKERRQ(ierr);
1194   ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),NULL);CHKERRQ(ierr);
1195   ierr = PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",mumps->id.CNTL(4),&mumps->id.CNTL(4),NULL);CHKERRQ(ierr);
1196   ierr = PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",mumps->id.CNTL(5),&mumps->id.CNTL(5),NULL);CHKERRQ(ierr);
1197   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);
1198 
1199   ierr = PetscOptionsString("-mat_mumps_ooc_tmpdir", "out of core directory", "None", mumps->id.ooc_tmpdir, mumps->id.ooc_tmpdir, 256, NULL);CHKERRQ(ierr);
1200 
1201   ierr = PetscOptionsIntArray("-mat_mumps_view_info","request INFO local to each processor","",info,&ninfo,NULL);CHKERRQ(ierr);
1202   if (ninfo) {
1203     if (ninfo > 40) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"number of INFO %d must <= 40\n",ninfo);
1204     ierr = PetscMalloc1(ninfo,&mumps->info);CHKERRQ(ierr);
1205     mumps->ninfo = ninfo;
1206     for (i=0; i<ninfo; i++) {
1207       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);
1208       else  mumps->info[i] = info[i];
1209     }
1210   }
1211 
1212   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1213   PetscFunctionReturn(0);
1214 }
1215 
1216 PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps)
1217 {
1218   PetscErrorCode ierr;
1219 
1220   PetscFunctionBegin;
1221   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)A), &mumps->myid);CHKERRQ(ierr);
1222   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&mumps->size);CHKERRQ(ierr);
1223   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(mumps->comm_mumps));CHKERRQ(ierr);
1224 
1225   mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps);
1226 
1227   mumps->id.job = JOB_INIT;
1228   mumps->id.par = 1;  /* host participates factorizaton and solve */
1229   mumps->id.sym = mumps->sym;
1230   PetscMUMPS_c(&mumps->id);
1231 
1232   mumps->scat_rhs     = NULL;
1233   mumps->scat_sol     = NULL;
1234 
1235   /* set PETSc-MUMPS default options - override MUMPS default */
1236   mumps->id.ICNTL(3) = 0;
1237   mumps->id.ICNTL(4) = 0;
1238   if (mumps->size == 1) {
1239     mumps->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
1240   } else {
1241     mumps->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
1242     mumps->id.ICNTL(20) = 0;   /* rhs is in dense format */
1243     mumps->id.ICNTL(21) = 1;   /* distributed solution */
1244   }
1245 
1246   /* schur */
1247   mumps->id.size_schur      = 0;
1248   mumps->id.listvar_schur   = NULL;
1249   mumps->id.schur           = NULL;
1250   mumps->sizeredrhs         = 0;
1251   mumps->schur_sol          = NULL;
1252   mumps->schur_sizesol      = 0;
1253   PetscFunctionReturn(0);
1254 }
1255 
1256 PetscErrorCode MatFactorSymbolic_MUMPS_ReportIfError(Mat F,Mat A,const MatFactorInfo *info,Mat_MUMPS *mumps)
1257 {
1258   PetscErrorCode ierr;
1259 
1260   PetscFunctionBegin;
1261   if (mumps->id.INFOG(1) < 0) {
1262     if (A->erroriffailure) {
1263       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
1264     } else {
1265       if (mumps->id.INFOG(1) == -6) {
1266         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);
1267         F->factorerrortype = MAT_FACTOR_STRUCT_ZEROPIVOT;
1268       } else if (mumps->id.INFOG(1) == -5 || mumps->id.INFOG(1) == -7) {
1269         ierr = PetscInfo2(F,"problem of workspace, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
1270         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1271       } else {
1272         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);
1273         F->factorerrortype = MAT_FACTOR_OTHER;
1274       }
1275     }
1276   }
1277   PetscFunctionReturn(0);
1278 }
1279 
1280 /* Note Petsc r(=c) permutation is used when mumps->id.ICNTL(7)==1 with centralized assembled matrix input; otherwise r and c are ignored */
1281 PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1282 {
1283   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1284   PetscErrorCode ierr;
1285   Vec            b;
1286   IS             is_iden;
1287   const PetscInt M = A->rmap->N;
1288 
1289   PetscFunctionBegin;
1290   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1291 
1292   /* Set MUMPS options from the options database */
1293   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1294 
1295   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1296 
1297   /* analysis phase */
1298   /*----------------*/
1299   mumps->id.job = JOB_FACTSYMBOLIC;
1300   mumps->id.n   = M;
1301   switch (mumps->id.ICNTL(18)) {
1302   case 0:  /* centralized assembled matrix input */
1303     if (!mumps->myid) {
1304       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1305       if (mumps->id.ICNTL(6)>1) {
1306         mumps->id.a = (MumpsScalar*)mumps->val;
1307       }
1308       if (mumps->id.ICNTL(7) == 1) { /* use user-provide matrix ordering - assuming r = c ordering */
1309         /*
1310         PetscBool      flag;
1311         ierr = ISEqual(r,c,&flag);CHKERRQ(ierr);
1312         if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"row_perm != col_perm");
1313         ierr = ISView(r,PETSC_VIEWER_STDOUT_SELF);
1314          */
1315         if (!mumps->myid) {
1316           const PetscInt *idx;
1317           PetscInt       i,*perm_in;
1318 
1319           ierr = PetscMalloc1(M,&perm_in);CHKERRQ(ierr);
1320           ierr = ISGetIndices(r,&idx);CHKERRQ(ierr);
1321 
1322           mumps->id.perm_in = perm_in;
1323           for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */
1324           ierr = ISRestoreIndices(r,&idx);CHKERRQ(ierr);
1325         }
1326       }
1327     }
1328     break;
1329   case 3:  /* distributed assembled matrix input (size>1) */
1330     mumps->id.nz_loc = mumps->nz;
1331     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1332     if (mumps->id.ICNTL(6)>1) {
1333       mumps->id.a_loc = (MumpsScalar*)mumps->val;
1334     }
1335     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1336     if (!mumps->myid) {
1337       ierr = VecCreateSeq(PETSC_COMM_SELF,A->rmap->N,&mumps->b_seq);CHKERRQ(ierr);
1338       ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->N,0,1,&is_iden);CHKERRQ(ierr);
1339     } else {
1340       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
1341       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
1342     }
1343     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1344     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
1345     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
1346     ierr = VecDestroy(&b);CHKERRQ(ierr);
1347     break;
1348   }
1349   PetscMUMPS_c(&mumps->id);
1350   ierr = MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);CHKERRQ(ierr);
1351 
1352   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1353   F->ops->solve           = MatSolve_MUMPS;
1354   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1355   F->ops->matsolve        = MatMatSolve_MUMPS;
1356   PetscFunctionReturn(0);
1357 }
1358 
1359 /* Note the Petsc r and c permutations are ignored */
1360 PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1361 {
1362   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1363   PetscErrorCode ierr;
1364   Vec            b;
1365   IS             is_iden;
1366   const PetscInt M = A->rmap->N;
1367 
1368   PetscFunctionBegin;
1369   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1370 
1371   /* Set MUMPS options from the options database */
1372   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1373 
1374   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1375 
1376   /* analysis phase */
1377   /*----------------*/
1378   mumps->id.job = JOB_FACTSYMBOLIC;
1379   mumps->id.n   = M;
1380   switch (mumps->id.ICNTL(18)) {
1381   case 0:  /* centralized assembled matrix input */
1382     if (!mumps->myid) {
1383       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1384       if (mumps->id.ICNTL(6)>1) {
1385         mumps->id.a = (MumpsScalar*)mumps->val;
1386       }
1387     }
1388     break;
1389   case 3:  /* distributed assembled matrix input (size>1) */
1390     mumps->id.nz_loc = mumps->nz;
1391     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1392     if (mumps->id.ICNTL(6)>1) {
1393       mumps->id.a_loc = (MumpsScalar*)mumps->val;
1394     }
1395     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1396     if (!mumps->myid) {
1397       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
1398       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
1399     } else {
1400       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
1401       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
1402     }
1403     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1404     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
1405     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
1406     ierr = VecDestroy(&b);CHKERRQ(ierr);
1407     break;
1408   }
1409   PetscMUMPS_c(&mumps->id);
1410   ierr = MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);CHKERRQ(ierr);
1411 
1412   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1413   F->ops->solve           = MatSolve_MUMPS;
1414   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1415   PetscFunctionReturn(0);
1416 }
1417 
1418 /* Note the Petsc r permutation and factor info are ignored */
1419 PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
1420 {
1421   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1422   PetscErrorCode ierr;
1423   Vec            b;
1424   IS             is_iden;
1425   const PetscInt M = A->rmap->N;
1426 
1427   PetscFunctionBegin;
1428   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1429 
1430   /* Set MUMPS options from the options database */
1431   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1432 
1433   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1434 
1435   /* analysis phase */
1436   /*----------------*/
1437   mumps->id.job = JOB_FACTSYMBOLIC;
1438   mumps->id.n   = M;
1439   switch (mumps->id.ICNTL(18)) {
1440   case 0:  /* centralized assembled matrix input */
1441     if (!mumps->myid) {
1442       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1443       if (mumps->id.ICNTL(6)>1) {
1444         mumps->id.a = (MumpsScalar*)mumps->val;
1445       }
1446     }
1447     break;
1448   case 3:  /* distributed assembled matrix input (size>1) */
1449     mumps->id.nz_loc = mumps->nz;
1450     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1451     if (mumps->id.ICNTL(6)>1) {
1452       mumps->id.a_loc = (MumpsScalar*)mumps->val;
1453     }
1454     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1455     if (!mumps->myid) {
1456       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
1457       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
1458     } else {
1459       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
1460       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
1461     }
1462     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1463     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
1464     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
1465     ierr = VecDestroy(&b);CHKERRQ(ierr);
1466     break;
1467   }
1468   PetscMUMPS_c(&mumps->id);
1469   ierr = MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);CHKERRQ(ierr);
1470 
1471   F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
1472   F->ops->solve                 = MatSolve_MUMPS;
1473   F->ops->solvetranspose        = MatSolve_MUMPS;
1474   F->ops->matsolve              = MatMatSolve_MUMPS;
1475 #if defined(PETSC_USE_COMPLEX)
1476   F->ops->getinertia = NULL;
1477 #else
1478   F->ops->getinertia = MatGetInertia_SBAIJMUMPS;
1479 #endif
1480   PetscFunctionReturn(0);
1481 }
1482 
1483 PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
1484 {
1485   PetscErrorCode    ierr;
1486   PetscBool         iascii;
1487   PetscViewerFormat format;
1488   Mat_MUMPS         *mumps=(Mat_MUMPS*)A->data;
1489 
1490   PetscFunctionBegin;
1491   /* check if matrix is mumps type */
1492   if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0);
1493 
1494   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1495   if (iascii) {
1496     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1497     if (format == PETSC_VIEWER_ASCII_INFO) {
1498       ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
1499       ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                   %d \n",mumps->id.sym);CHKERRQ(ierr);
1500       ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):            %d \n",mumps->id.par);CHKERRQ(ierr);
1501       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):         %d \n",mumps->id.ICNTL(1));CHKERRQ(ierr);
1502       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg): %d \n",mumps->id.ICNTL(2));CHKERRQ(ierr);
1503       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):   %d \n",mumps->id.ICNTL(3));CHKERRQ(ierr);
1504       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):        %d \n",mumps->id.ICNTL(4));CHKERRQ(ierr);
1505       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):         %d \n",mumps->id.ICNTL(5));CHKERRQ(ierr);
1506       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):        %d \n",mumps->id.ICNTL(6));CHKERRQ(ierr);
1507       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (sequential matrix ordering):%d \n",mumps->id.ICNTL(7));CHKERRQ(ierr);
1508       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scaling strategy):        %d \n",mumps->id.ICNTL(8));CHKERRQ(ierr);
1509       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements):  %d \n",mumps->id.ICNTL(10));CHKERRQ(ierr);
1510       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):          %d \n",mumps->id.ICNTL(11));CHKERRQ(ierr);
1511       if (mumps->id.ICNTL(11)>0) {
1512         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(4) (inf norm of input mat):        %g\n",mumps->id.RINFOG(4));CHKERRQ(ierr);
1513         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(5) (inf norm of solution):         %g\n",mumps->id.RINFOG(5));CHKERRQ(ierr);
1514         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(6) (inf norm of residual):         %g\n",mumps->id.RINFOG(6));CHKERRQ(ierr);
1515         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8));CHKERRQ(ierr);
1516         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(9) (error estimate):               %g \n",mumps->id.RINFOG(9));CHKERRQ(ierr);
1517         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11));CHKERRQ(ierr);
1518       }
1519       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",mumps->id.ICNTL(12));CHKERRQ(ierr);
1520       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",mumps->id.ICNTL(13));CHKERRQ(ierr);
1521       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",mumps->id.ICNTL(14));CHKERRQ(ierr);
1522       /* ICNTL(15-17) not used */
1523       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",mumps->id.ICNTL(18));CHKERRQ(ierr);
1524       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Schur complement info):                       %d \n",mumps->id.ICNTL(19));CHKERRQ(ierr);
1525       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",mumps->id.ICNTL(20));CHKERRQ(ierr);
1526       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d \n",mumps->id.ICNTL(21));CHKERRQ(ierr);
1527       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",mumps->id.ICNTL(22));CHKERRQ(ierr);
1528       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",mumps->id.ICNTL(23));CHKERRQ(ierr);
1529 
1530       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",mumps->id.ICNTL(24));CHKERRQ(ierr);
1531       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",mumps->id.ICNTL(25));CHKERRQ(ierr);
1532       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",mumps->id.ICNTL(26));CHKERRQ(ierr);
1533       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",mumps->id.ICNTL(27));CHKERRQ(ierr);
1534       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(28) (use parallel or sequential ordering):        %d \n",mumps->id.ICNTL(28));CHKERRQ(ierr);
1535       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(29) (parallel ordering):                          %d \n",mumps->id.ICNTL(29));CHKERRQ(ierr);
1536 
1537       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(30) (user-specified set of entries in inv(A)):    %d \n",mumps->id.ICNTL(30));CHKERRQ(ierr);
1538       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(31) (factors is discarded in the solve phase):    %d \n",mumps->id.ICNTL(31));CHKERRQ(ierr);
1539       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(33) (compute determinant):                        %d \n",mumps->id.ICNTL(33));CHKERRQ(ierr);
1540       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(35) (activate BLR based factorization):           %d \n",mumps->id.ICNTL(35));CHKERRQ(ierr);
1541 
1542       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",mumps->id.CNTL(1));CHKERRQ(ierr);
1543       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2));CHKERRQ(ierr);
1544       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",mumps->id.CNTL(3));CHKERRQ(ierr);
1545       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",mumps->id.CNTL(4));CHKERRQ(ierr);
1546       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",mumps->id.CNTL(5));CHKERRQ(ierr);
1547       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(7) (dropping parameter for BLR):       %g \n",mumps->id.CNTL(7));CHKERRQ(ierr);
1548 
1549       /* infomation local to each processor */
1550       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);
1551       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
1552       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %g \n",mumps->myid,mumps->id.RINFO(1));CHKERRQ(ierr);
1553       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1554       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);
1555       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(2));CHKERRQ(ierr);
1556       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1557       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);
1558       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(3));CHKERRQ(ierr);
1559       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1560 
1561       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr);
1562       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] %d \n",mumps->myid,mumps->id.INFO(15));CHKERRQ(ierr);
1563       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1564 
1565       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr);
1566       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(16));CHKERRQ(ierr);
1567       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1568 
1569       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);
1570       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(23));CHKERRQ(ierr);
1571       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1572 
1573       if (mumps->ninfo && mumps->ninfo <= 40){
1574         PetscInt i;
1575         for (i=0; i<mumps->ninfo; i++){
1576           ierr = PetscViewerASCIIPrintf(viewer, "  INFO(%d): \n",mumps->info[i]);CHKERRQ(ierr);
1577           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(mumps->info[i]));CHKERRQ(ierr);
1578           ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1579         }
1580       }
1581 
1582 
1583       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
1584 
1585       if (!mumps->myid) { /* information from the host */
1586         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",mumps->id.RINFOG(1));CHKERRQ(ierr);
1587         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",mumps->id.RINFOG(2));CHKERRQ(ierr);
1588         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",mumps->id.RINFOG(3));CHKERRQ(ierr);
1589         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);
1590 
1591         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(3));CHKERRQ(ierr);
1592         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(4));CHKERRQ(ierr);
1593         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",mumps->id.INFOG(5));CHKERRQ(ierr);
1594         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",mumps->id.INFOG(6));CHKERRQ(ierr);
1595         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively use after analysis): %d \n",mumps->id.INFOG(7));CHKERRQ(ierr);
1596         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",mumps->id.INFOG(8));CHKERRQ(ierr);
1597         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",mumps->id.INFOG(9));CHKERRQ(ierr);
1598         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",mumps->id.INFOG(10));CHKERRQ(ierr);
1599         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",mumps->id.INFOG(11));CHKERRQ(ierr);
1600         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",mumps->id.INFOG(12));CHKERRQ(ierr);
1601         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",mumps->id.INFOG(13));CHKERRQ(ierr);
1602         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",mumps->id.INFOG(14));CHKERRQ(ierr);
1603         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",mumps->id.INFOG(15));CHKERRQ(ierr);
1604         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);
1605         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);
1606         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);
1607         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);
1608         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",mumps->id.INFOG(20));CHKERRQ(ierr);
1609         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);
1610         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);
1611         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",mumps->id.INFOG(23));CHKERRQ(ierr);
1612         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",mumps->id.INFOG(24));CHKERRQ(ierr);
1613         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",mumps->id.INFOG(25));CHKERRQ(ierr);
1614         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(28) (after factorization: number of null pivots encountered): %d\n",mumps->id.INFOG(28));CHKERRQ(ierr);
1615         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);
1616         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);
1617         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(32) (after analysis: type of analysis done): %d\n",mumps->id.INFOG(32));CHKERRQ(ierr);
1618         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(33) (value used for ICNTL(8)): %d\n",mumps->id.INFOG(33));CHKERRQ(ierr);
1619         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(34) (exponent of the determinant if determinant is requested): %d\n",mumps->id.INFOG(34));CHKERRQ(ierr);
1620       }
1621     }
1622   }
1623   PetscFunctionReturn(0);
1624 }
1625 
1626 PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
1627 {
1628   Mat_MUMPS *mumps =(Mat_MUMPS*)A->data;
1629 
1630   PetscFunctionBegin;
1631   info->block_size        = 1.0;
1632   info->nz_allocated      = mumps->id.INFOG(20);
1633   info->nz_used           = mumps->id.INFOG(20);
1634   info->nz_unneeded       = 0.0;
1635   info->assemblies        = 0.0;
1636   info->mallocs           = 0.0;
1637   info->memory            = 0.0;
1638   info->fill_ratio_given  = 0;
1639   info->fill_ratio_needed = 0;
1640   info->factor_mallocs    = 0;
1641   PetscFunctionReturn(0);
1642 }
1643 
1644 /* -------------------------------------------------------------------------------------------*/
1645 PetscErrorCode MatFactorSetSchurIS_MUMPS(Mat F, IS is)
1646 {
1647   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1648   const PetscInt *idxs;
1649   PetscInt       size,i;
1650   PetscErrorCode ierr;
1651 
1652   PetscFunctionBegin;
1653   ierr = ISGetLocalSize(is,&size);CHKERRQ(ierr);
1654   if (mumps->size > 1) {
1655     PetscBool ls,gs;
1656 
1657     ls   = mumps->myid ? (size ? PETSC_FALSE : PETSC_TRUE) : PETSC_TRUE;
1658     ierr = MPI_Allreduce(&ls,&gs,1,MPIU_BOOL,MPI_LAND,mumps->comm_mumps);CHKERRQ(ierr);
1659     if (!gs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MUMPS distributed parallel Schur complements not yet supported from PETSc\n");
1660   }
1661   if (mumps->id.size_schur != size) {
1662     ierr = PetscFree2(mumps->id.listvar_schur,mumps->id.schur);CHKERRQ(ierr);
1663     mumps->id.size_schur = size;
1664     mumps->id.schur_lld  = size;
1665     ierr = PetscMalloc2(size,&mumps->id.listvar_schur,size*size,&mumps->id.schur);CHKERRQ(ierr);
1666   }
1667 
1668   /* Schur complement matrix */
1669   ierr = MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.size_schur,(PetscScalar*)mumps->id.schur,&F->schur);CHKERRQ(ierr);
1670   if (mumps->sym == 1) {
1671     ierr = MatSetOption(F->schur,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr);
1672   }
1673 
1674   /* MUMPS expects Fortran style indices */
1675   ierr = ISGetIndices(is,&idxs);CHKERRQ(ierr);
1676   ierr = PetscMemcpy(mumps->id.listvar_schur,idxs,size*sizeof(PetscInt));CHKERRQ(ierr);
1677   for (i=0;i<size;i++) mumps->id.listvar_schur[i]++;
1678   ierr = ISRestoreIndices(is,&idxs);CHKERRQ(ierr);
1679   if (mumps->size > 1) {
1680     mumps->id.ICNTL(19) = 1; /* MUMPS returns Schur centralized on the host */
1681   } else {
1682     if (F->factortype == MAT_FACTOR_LU) {
1683       mumps->id.ICNTL(19) = 3; /* MUMPS returns full matrix */
1684     } else {
1685       mumps->id.ICNTL(19) = 2; /* MUMPS returns lower triangular part */
1686     }
1687   }
1688   /* set a special value of ICNTL (not handled my MUMPS) to be used in the solve phase by PETSc */
1689   mumps->id.ICNTL(26) = -1;
1690   PetscFunctionReturn(0);
1691 }
1692 
1693 /* -------------------------------------------------------------------------------------------*/
1694 PetscErrorCode MatFactorCreateSchurComplement_MUMPS(Mat F,Mat* S)
1695 {
1696   Mat            St;
1697   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1698   PetscScalar    *array;
1699 #if defined(PETSC_USE_COMPLEX)
1700   PetscScalar    im = PetscSqrtScalar((PetscScalar)-1.0);
1701 #endif
1702   PetscErrorCode ierr;
1703 
1704   PetscFunctionBegin;
1705   if (!mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it");
1706   ierr = MatCreate(PETSC_COMM_SELF,&St);CHKERRQ(ierr);
1707   ierr = MatSetSizes(St,PETSC_DECIDE,PETSC_DECIDE,mumps->id.size_schur,mumps->id.size_schur);CHKERRQ(ierr);
1708   ierr = MatSetType(St,MATDENSE);CHKERRQ(ierr);
1709   ierr = MatSetUp(St);CHKERRQ(ierr);
1710   ierr = MatDenseGetArray(St,&array);CHKERRQ(ierr);
1711   if (!mumps->sym) { /* MUMPS always return a full matrix */
1712     if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
1713       PetscInt i,j,N=mumps->id.size_schur;
1714       for (i=0;i<N;i++) {
1715         for (j=0;j<N;j++) {
1716 #if !defined(PETSC_USE_COMPLEX)
1717           PetscScalar val = mumps->id.schur[i*N+j];
1718 #else
1719           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
1720 #endif
1721           array[j*N+i] = val;
1722         }
1723       }
1724     } else { /* stored by columns */
1725       ierr = PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
1726     }
1727   } else { /* either full or lower-triangular (not packed) */
1728     if (mumps->id.ICNTL(19) == 2) { /* lower triangular stored by columns */
1729       PetscInt i,j,N=mumps->id.size_schur;
1730       for (i=0;i<N;i++) {
1731         for (j=i;j<N;j++) {
1732 #if !defined(PETSC_USE_COMPLEX)
1733           PetscScalar val = mumps->id.schur[i*N+j];
1734 #else
1735           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
1736 #endif
1737           array[i*N+j] = val;
1738           array[j*N+i] = val;
1739         }
1740       }
1741     } else if (mumps->id.ICNTL(19) == 3) { /* full matrix */
1742       ierr = PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
1743     } else { /* ICNTL(19) == 1 lower triangular stored by rows */
1744       PetscInt i,j,N=mumps->id.size_schur;
1745       for (i=0;i<N;i++) {
1746         for (j=0;j<i+1;j++) {
1747 #if !defined(PETSC_USE_COMPLEX)
1748           PetscScalar val = mumps->id.schur[i*N+j];
1749 #else
1750           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
1751 #endif
1752           array[i*N+j] = val;
1753           array[j*N+i] = val;
1754         }
1755       }
1756     }
1757   }
1758   ierr = MatDenseRestoreArray(St,&array);CHKERRQ(ierr);
1759   *S   = St;
1760   PetscFunctionReturn(0);
1761 }
1762 
1763 /* -------------------------------------------------------------------------------------------*/
1764 PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival)
1765 {
1766   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
1767 
1768   PetscFunctionBegin;
1769   mumps->id.ICNTL(icntl) = ival;
1770   PetscFunctionReturn(0);
1771 }
1772 
1773 PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt *ival)
1774 {
1775   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
1776 
1777   PetscFunctionBegin;
1778   *ival = mumps->id.ICNTL(icntl);
1779   PetscFunctionReturn(0);
1780 }
1781 
1782 /*@
1783   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()
1784 
1785    Logically Collective on Mat
1786 
1787    Input Parameters:
1788 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1789 .  icntl - index of MUMPS parameter array ICNTL()
1790 -  ival - value of MUMPS ICNTL(icntl)
1791 
1792   Options Database:
1793 .   -mat_mumps_icntl_<icntl> <ival>
1794 
1795    Level: beginner
1796 
1797    References:
1798 .     MUMPS Users' Guide
1799 
1800 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
1801  @*/
1802 PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival)
1803 {
1804   PetscErrorCode ierr;
1805 
1806   PetscFunctionBegin;
1807   PetscValidType(F,1);
1808   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
1809   PetscValidLogicalCollectiveInt(F,icntl,2);
1810   PetscValidLogicalCollectiveInt(F,ival,3);
1811   ierr = PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr);
1812   PetscFunctionReturn(0);
1813 }
1814 
1815 /*@
1816   MatMumpsGetIcntl - Get MUMPS parameter ICNTL()
1817 
1818    Logically Collective on Mat
1819 
1820    Input Parameters:
1821 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1822 -  icntl - index of MUMPS parameter array ICNTL()
1823 
1824   Output Parameter:
1825 .  ival - value of MUMPS ICNTL(icntl)
1826 
1827    Level: beginner
1828 
1829    References:
1830 .     MUMPS Users' Guide
1831 
1832 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
1833 @*/
1834 PetscErrorCode MatMumpsGetIcntl(Mat F,PetscInt icntl,PetscInt *ival)
1835 {
1836   PetscErrorCode ierr;
1837 
1838   PetscFunctionBegin;
1839   PetscValidType(F,1);
1840   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
1841   PetscValidLogicalCollectiveInt(F,icntl,2);
1842   PetscValidIntPointer(ival,3);
1843   ierr = PetscUseMethod(F,"MatMumpsGetIcntl_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
1844   PetscFunctionReturn(0);
1845 }
1846 
1847 /* -------------------------------------------------------------------------------------------*/
1848 PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal val)
1849 {
1850   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
1851 
1852   PetscFunctionBegin;
1853   mumps->id.CNTL(icntl) = val;
1854   PetscFunctionReturn(0);
1855 }
1856 
1857 PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal *val)
1858 {
1859   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
1860 
1861   PetscFunctionBegin;
1862   *val = mumps->id.CNTL(icntl);
1863   PetscFunctionReturn(0);
1864 }
1865 
1866 /*@
1867   MatMumpsSetCntl - Set MUMPS parameter CNTL()
1868 
1869    Logically Collective on Mat
1870 
1871    Input Parameters:
1872 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1873 .  icntl - index of MUMPS parameter array CNTL()
1874 -  val - value of MUMPS CNTL(icntl)
1875 
1876   Options Database:
1877 .   -mat_mumps_cntl_<icntl> <val>
1878 
1879    Level: beginner
1880 
1881    References:
1882 .     MUMPS Users' Guide
1883 
1884 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
1885 @*/
1886 PetscErrorCode MatMumpsSetCntl(Mat F,PetscInt icntl,PetscReal val)
1887 {
1888   PetscErrorCode ierr;
1889 
1890   PetscFunctionBegin;
1891   PetscValidType(F,1);
1892   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
1893   PetscValidLogicalCollectiveInt(F,icntl,2);
1894   PetscValidLogicalCollectiveReal(F,val,3);
1895   ierr = PetscTryMethod(F,"MatMumpsSetCntl_C",(Mat,PetscInt,PetscReal),(F,icntl,val));CHKERRQ(ierr);
1896   PetscFunctionReturn(0);
1897 }
1898 
1899 /*@
1900   MatMumpsGetCntl - Get MUMPS parameter CNTL()
1901 
1902    Logically Collective on Mat
1903 
1904    Input Parameters:
1905 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1906 -  icntl - index of MUMPS parameter array CNTL()
1907 
1908   Output Parameter:
1909 .  val - value of MUMPS CNTL(icntl)
1910 
1911    Level: beginner
1912 
1913    References:
1914 .      MUMPS Users' Guide
1915 
1916 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
1917 @*/
1918 PetscErrorCode MatMumpsGetCntl(Mat F,PetscInt icntl,PetscReal *val)
1919 {
1920   PetscErrorCode ierr;
1921 
1922   PetscFunctionBegin;
1923   PetscValidType(F,1);
1924   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
1925   PetscValidLogicalCollectiveInt(F,icntl,2);
1926   PetscValidRealPointer(val,3);
1927   ierr = PetscUseMethod(F,"MatMumpsGetCntl_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
1928   PetscFunctionReturn(0);
1929 }
1930 
1931 PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F,PetscInt icntl,PetscInt *info)
1932 {
1933   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
1934 
1935   PetscFunctionBegin;
1936   *info = mumps->id.INFO(icntl);
1937   PetscFunctionReturn(0);
1938 }
1939 
1940 PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F,PetscInt icntl,PetscInt *infog)
1941 {
1942   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
1943 
1944   PetscFunctionBegin;
1945   *infog = mumps->id.INFOG(icntl);
1946   PetscFunctionReturn(0);
1947 }
1948 
1949 PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfo)
1950 {
1951   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
1952 
1953   PetscFunctionBegin;
1954   *rinfo = mumps->id.RINFO(icntl);
1955   PetscFunctionReturn(0);
1956 }
1957 
1958 PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfog)
1959 {
1960   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
1961 
1962   PetscFunctionBegin;
1963   *rinfog = mumps->id.RINFOG(icntl);
1964   PetscFunctionReturn(0);
1965 }
1966 
1967 /*@
1968   MatMumpsGetInfo - Get MUMPS parameter INFO()
1969 
1970    Logically Collective on Mat
1971 
1972    Input Parameters:
1973 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1974 -  icntl - index of MUMPS parameter array INFO()
1975 
1976   Output Parameter:
1977 .  ival - value of MUMPS INFO(icntl)
1978 
1979    Level: beginner
1980 
1981    References:
1982 .      MUMPS Users' Guide
1983 
1984 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
1985 @*/
1986 PetscErrorCode MatMumpsGetInfo(Mat F,PetscInt icntl,PetscInt *ival)
1987 {
1988   PetscErrorCode ierr;
1989 
1990   PetscFunctionBegin;
1991   PetscValidType(F,1);
1992   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
1993   PetscValidIntPointer(ival,3);
1994   ierr = PetscUseMethod(F,"MatMumpsGetInfo_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
1995   PetscFunctionReturn(0);
1996 }
1997 
1998 /*@
1999   MatMumpsGetInfog - Get MUMPS parameter INFOG()
2000 
2001    Logically Collective on Mat
2002 
2003    Input Parameters:
2004 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2005 -  icntl - index of MUMPS parameter array INFOG()
2006 
2007   Output Parameter:
2008 .  ival - value of MUMPS INFOG(icntl)
2009 
2010    Level: beginner
2011 
2012    References:
2013 .      MUMPS Users' Guide
2014 
2015 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2016 @*/
2017 PetscErrorCode MatMumpsGetInfog(Mat F,PetscInt icntl,PetscInt *ival)
2018 {
2019   PetscErrorCode ierr;
2020 
2021   PetscFunctionBegin;
2022   PetscValidType(F,1);
2023   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2024   PetscValidIntPointer(ival,3);
2025   ierr = PetscUseMethod(F,"MatMumpsGetInfog_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
2026   PetscFunctionReturn(0);
2027 }
2028 
2029 /*@
2030   MatMumpsGetRinfo - Get MUMPS parameter RINFO()
2031 
2032    Logically Collective on Mat
2033 
2034    Input Parameters:
2035 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2036 -  icntl - index of MUMPS parameter array RINFO()
2037 
2038   Output Parameter:
2039 .  val - value of MUMPS RINFO(icntl)
2040 
2041    Level: beginner
2042 
2043    References:
2044 .       MUMPS Users' Guide
2045 
2046 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2047 @*/
2048 PetscErrorCode MatMumpsGetRinfo(Mat F,PetscInt icntl,PetscReal *val)
2049 {
2050   PetscErrorCode ierr;
2051 
2052   PetscFunctionBegin;
2053   PetscValidType(F,1);
2054   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2055   PetscValidRealPointer(val,3);
2056   ierr = PetscUseMethod(F,"MatMumpsGetRinfo_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
2057   PetscFunctionReturn(0);
2058 }
2059 
2060 /*@
2061   MatMumpsGetRinfog - Get MUMPS parameter RINFOG()
2062 
2063    Logically Collective on Mat
2064 
2065    Input Parameters:
2066 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2067 -  icntl - index of MUMPS parameter array RINFOG()
2068 
2069   Output Parameter:
2070 .  val - value of MUMPS RINFOG(icntl)
2071 
2072    Level: beginner
2073 
2074    References:
2075 .      MUMPS Users' Guide
2076 
2077 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2078 @*/
2079 PetscErrorCode MatMumpsGetRinfog(Mat F,PetscInt icntl,PetscReal *val)
2080 {
2081   PetscErrorCode ierr;
2082 
2083   PetscFunctionBegin;
2084   PetscValidType(F,1);
2085   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2086   PetscValidRealPointer(val,3);
2087   ierr = PetscUseMethod(F,"MatMumpsGetRinfog_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
2088   PetscFunctionReturn(0);
2089 }
2090 
2091 /*MC
2092   MATSOLVERMUMPS -  A matrix type providing direct solvers (LU and Cholesky) for
2093   distributed and sequential matrices via the external package MUMPS.
2094 
2095   Works with MATAIJ and MATSBAIJ matrices
2096 
2097   Use ./configure --download-mumps --download-scalapack --download-parmetis --download-metis --download-ptscotch  to have PETSc installed with MUMPS
2098 
2099   Use -pc_type cholesky or lu -pc_factor_mat_solver_type mumps to use this direct solver
2100 
2101   Options Database Keys:
2102 +  -mat_mumps_icntl_1 - ICNTL(1): output stream for error messages
2103 .  -mat_mumps_icntl_2 - ICNTL(2): output stream for diagnostic printing, statistics, and warning
2104 .  -mat_mumps_icntl_3 -  ICNTL(3): output stream for global information, collected on the host
2105 .  -mat_mumps_icntl_4 -  ICNTL(4): level of printing (0 to 4)
2106 .  -mat_mumps_icntl_6 - ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7)
2107 .  -mat_mumps_icntl_7 - ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis
2108 .  -mat_mumps_icntl_8  - ICNTL(8): scaling strategy (-2 to 8 or 77)
2109 .  -mat_mumps_icntl_10  - ICNTL(10): max num of refinements
2110 .  -mat_mumps_icntl_11  - ICNTL(11): statistics related to an error analysis (via -ksp_view)
2111 .  -mat_mumps_icntl_12  - ICNTL(12): an ordering strategy for symmetric matrices (0 to 3)
2112 .  -mat_mumps_icntl_13  - ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting
2113 .  -mat_mumps_icntl_14  - ICNTL(14): percentage increase in the estimated working space
2114 .  -mat_mumps_icntl_19  - ICNTL(19): computes the Schur complement
2115 .  -mat_mumps_icntl_22  - ICNTL(22): in-core/out-of-core factorization and solve (0 or 1)
2116 .  -mat_mumps_icntl_23  - ICNTL(23): max size of the working memory (MB) that can allocate per processor
2117 .  -mat_mumps_icntl_24  - ICNTL(24): detection of null pivot rows (0 or 1)
2118 .  -mat_mumps_icntl_25  - ICNTL(25): compute a solution of a deficient matrix and a null space basis
2119 .  -mat_mumps_icntl_26  - ICNTL(26): drives the solution phase if a Schur complement matrix
2120 .  -mat_mumps_icntl_28  - ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering
2121 .  -mat_mumps_icntl_29 - ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis
2122 .  -mat_mumps_icntl_30 - ICNTL(30): compute user-specified set of entries in inv(A)
2123 .  -mat_mumps_icntl_31 - ICNTL(31): indicates which factors may be discarded during factorization
2124 .  -mat_mumps_icntl_33 - ICNTL(33): compute determinant
2125 .  -mat_mumps_cntl_1  - CNTL(1): relative pivoting threshold
2126 .  -mat_mumps_cntl_2  -  CNTL(2): stopping criterion of refinement
2127 .  -mat_mumps_cntl_3 - CNTL(3): absolute pivoting threshold
2128 .  -mat_mumps_cntl_4 - CNTL(4): value for static pivoting
2129 -  -mat_mumps_cntl_5 - CNTL(5): fixation for null pivots
2130 
2131   Level: beginner
2132 
2133     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
2134 $          KSPGetPC(ksp,&pc);
2135 $          PCFactorGetMatrix(pc,&mat);
2136 $          MatMumpsGetInfo(mat,....);
2137 $          MatMumpsGetInfog(mat,....); etc.
2138            Or you can run with -ksp_error_if_not_converged and the program will be stopped and the information printed in the error message.
2139 
2140 .seealso: PCFactorSetMatSolverType(), MatSolverType, MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog(), KSPGetPC(), PCGetFactor(), PCFactorGetMatrix()
2141 
2142 M*/
2143 
2144 static PetscErrorCode MatFactorGetSolverType_mumps(Mat A,MatSolverType *type)
2145 {
2146   PetscFunctionBegin;
2147   *type = MATSOLVERMUMPS;
2148   PetscFunctionReturn(0);
2149 }
2150 
2151 /* MatGetFactor for Seq and MPI AIJ matrices */
2152 static PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
2153 {
2154   Mat            B;
2155   PetscErrorCode ierr;
2156   Mat_MUMPS      *mumps;
2157   PetscBool      isSeqAIJ;
2158 
2159   PetscFunctionBegin;
2160   /* Create the factorization matrix */
2161   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr);
2162   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
2163   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
2164   ierr = PetscStrallocpy("mumps",&((PetscObject)B)->type_name);CHKERRQ(ierr);
2165   ierr = MatSetUp(B);CHKERRQ(ierr);
2166 
2167   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
2168 
2169   B->ops->view        = MatView_MUMPS;
2170   B->ops->getinfo     = MatGetInfo_MUMPS;
2171 
2172   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);CHKERRQ(ierr);
2173   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr);
2174   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr);
2175   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
2176   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
2177   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
2178   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
2179   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
2180   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
2181   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
2182   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
2183 
2184   if (ftype == MAT_FACTOR_LU) {
2185     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
2186     B->factortype            = MAT_FACTOR_LU;
2187     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
2188     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
2189     mumps->sym = 0;
2190   } else {
2191     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2192     B->factortype                  = MAT_FACTOR_CHOLESKY;
2193     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
2194     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
2195 #if defined(PETSC_USE_COMPLEX)
2196     mumps->sym = 2;
2197 #else
2198     if (A->spd_set && A->spd) mumps->sym = 1;
2199     else                      mumps->sym = 2;
2200 #endif
2201   }
2202 
2203   /* set solvertype */
2204   ierr = PetscFree(B->solvertype);CHKERRQ(ierr);
2205   ierr = PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);CHKERRQ(ierr);
2206 
2207   mumps->isAIJ    = PETSC_TRUE;
2208   B->ops->destroy = MatDestroy_MUMPS;
2209   B->data        = (void*)mumps;
2210 
2211   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
2212 
2213   *F = B;
2214   PetscFunctionReturn(0);
2215 }
2216 
2217 /* MatGetFactor for Seq and MPI SBAIJ matrices */
2218 static PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
2219 {
2220   Mat            B;
2221   PetscErrorCode ierr;
2222   Mat_MUMPS      *mumps;
2223   PetscBool      isSeqSBAIJ;
2224 
2225   PetscFunctionBegin;
2226   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
2227   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");
2228   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
2229   /* Create the factorization matrix */
2230   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
2231   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
2232   ierr = PetscStrallocpy("mumps",&((PetscObject)B)->type_name);CHKERRQ(ierr);
2233   ierr = MatSetUp(B);CHKERRQ(ierr);
2234 
2235   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
2236   if (isSeqSBAIJ) {
2237     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
2238   } else {
2239     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
2240   }
2241 
2242   B->ops->getinfo                = MatGetInfo_External;
2243   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2244   B->ops->view                   = MatView_MUMPS;
2245 
2246   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);CHKERRQ(ierr);
2247   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr);
2248   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr);
2249   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
2250   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
2251   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
2252   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
2253   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
2254   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
2255   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
2256   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
2257 
2258   B->factortype = MAT_FACTOR_CHOLESKY;
2259 #if defined(PETSC_USE_COMPLEX)
2260   mumps->sym = 2;
2261 #else
2262   if (A->spd_set && A->spd) mumps->sym = 1;
2263   else                      mumps->sym = 2;
2264 #endif
2265 
2266   /* set solvertype */
2267   ierr = PetscFree(B->solvertype);CHKERRQ(ierr);
2268   ierr = PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);CHKERRQ(ierr);
2269 
2270   mumps->isAIJ    = PETSC_FALSE;
2271   B->ops->destroy = MatDestroy_MUMPS;
2272   B->data        = (void*)mumps;
2273 
2274   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
2275 
2276   *F = B;
2277   PetscFunctionReturn(0);
2278 }
2279 
2280 static PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
2281 {
2282   Mat            B;
2283   PetscErrorCode ierr;
2284   Mat_MUMPS      *mumps;
2285   PetscBool      isSeqBAIJ;
2286 
2287   PetscFunctionBegin;
2288   /* Create the factorization matrix */
2289   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr);
2290   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
2291   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
2292   ierr = PetscStrallocpy("mumps",&((PetscObject)B)->type_name);CHKERRQ(ierr);
2293   ierr = MatSetUp(B);CHKERRQ(ierr);
2294 
2295   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
2296   if (ftype == MAT_FACTOR_LU) {
2297     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
2298     B->factortype            = MAT_FACTOR_LU;
2299     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
2300     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
2301     mumps->sym = 0;
2302   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n");
2303 
2304   B->ops->getinfo     = MatGetInfo_External;
2305   B->ops->view        = MatView_MUMPS;
2306 
2307   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);CHKERRQ(ierr);
2308   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr);
2309   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr);
2310   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
2311   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
2312   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
2313   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
2314   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
2315   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
2316   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
2317   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
2318 
2319   /* set solvertype */
2320   ierr = PetscFree(B->solvertype);CHKERRQ(ierr);
2321   ierr = PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);CHKERRQ(ierr);
2322 
2323   mumps->isAIJ    = PETSC_TRUE;
2324   B->ops->destroy = MatDestroy_MUMPS;
2325   B->data        = (void*)mumps;
2326 
2327   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
2328 
2329   *F = B;
2330   PetscFunctionReturn(0);
2331 }
2332 
2333 PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_MUMPS(void)
2334 {
2335   PetscErrorCode ierr;
2336 
2337   PetscFunctionBegin;
2338   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr);
2339   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr);
2340   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr);
2341   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr);
2342   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATMPISBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr);
2343   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr);
2344   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr);
2345   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr);
2346   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr);
2347   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr);
2348   PetscFunctionReturn(0);
2349 }
2350 
2351