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