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