xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision ca4a67af8ca7cf2d1ac25cfe8f5aed45e6b9324f)
1 #define PETSCMAT_DLL
2 
3 /*
4     Provides an interface to the MUMPS sparse solver
5 */
6 #include "../src/mat/impls/aij/seq/aij.h"
7 #include "../src/mat/impls/aij/mpi/mpiaij.h"
8 #include "../src/mat/impls/sbaij/seq/sbaij.h"
9 #include "../src/mat/impls/sbaij/mpi/mpisbaij.h"
10 
11 EXTERN_C_BEGIN
12 #if defined(PETSC_USE_COMPLEX)
13 #include "zmumps_c.h"
14 #else
15 #include "dmumps_c.h"
16 #endif
17 EXTERN_C_END
18 #define JOB_INIT -1
19 #define JOB_END -2
20 /* macros s.t. indices match MUMPS documentation */
21 #define ICNTL(I) icntl[(I)-1]
22 #define CNTL(I) cntl[(I)-1]
23 #define INFOG(I) infog[(I)-1]
24 #define INFO(I) info[(I)-1]
25 #define RINFOG(I) rinfog[(I)-1]
26 #define RINFO(I) rinfo[(I)-1]
27 
28 typedef struct {
29 #if defined(PETSC_USE_COMPLEX)
30   ZMUMPS_STRUC_C id;
31 #else
32   DMUMPS_STRUC_C id;
33 #endif
34   MatStructure   matstruc;
35   PetscMPIInt    myid,size;
36   PetscInt       *irn,*jcn,sym,nSolve;
37   PetscScalar    *val;
38   MPI_Comm       comm_mumps;
39   VecScatter     scat_rhs, scat_sol;
40   PetscTruth     isAIJ,CleanUpMUMPS;
41   Vec            b_seq,x_seq;
42   PetscErrorCode (*MatDestroy)(Mat);
43 } Mat_MUMPS;
44 
45 EXTERN PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);
46 
47 /* convert Petsc mpiaij matrix to triples: row[nz], col[nz], val[nz] */
48 /*
49   input:
50     A       - matrix in mpiaij or mpisbaij (bs=1) format
51     shift   - 0: C style output triple; 1: Fortran style output triple.
52     valOnly - FALSE: spaces are allocated and values are set for the triple
53               TRUE:  only the values in v array are updated
54   output:
55     nnz     - dim of r, c, and v (number of local nonzero entries of A)
56     r, c, v - row and col index, matrix values (matrix triples)
57  */
58 PetscErrorCode MatConvertToTriples(Mat A,int shift,PetscTruth valOnly,int *nnz,int **r, int **c, PetscScalar **v)
59 {
60   PetscInt       *ai, *aj, *bi, *bj, rstart,nz, *garray;
61   PetscErrorCode ierr;
62   PetscInt       i,j,jj,jB,irow,m=A->rmap->n,*ajj,*bjj,countA,countB,colA_start,jcol;
63   PetscInt       *row,*col;
64   PetscScalar    *av, *bv,*val;
65   PetscTruth     isAIJ;
66 
67   PetscFunctionBegin;
68   ierr = PetscTypeCompare((PetscObject)A,MATMPIAIJ,&isAIJ);CHKERRQ(ierr);
69   if (isAIJ){
70     Mat_MPIAIJ    *mat =  (Mat_MPIAIJ*)A->data;
71     Mat_SeqAIJ    *aa=(Mat_SeqAIJ*)(mat->A)->data;
72     Mat_SeqAIJ    *bb=(Mat_SeqAIJ*)(mat->B)->data;
73     nz = aa->nz + bb->nz;
74     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
75     garray = mat->garray;
76     av=aa->a; bv=bb->a;
77 
78   } else {
79     Mat_MPISBAIJ  *mat =  (Mat_MPISBAIJ*)A->data;
80     Mat_SeqSBAIJ  *aa=(Mat_SeqSBAIJ*)(mat->A)->data;
81     Mat_SeqBAIJ    *bb=(Mat_SeqBAIJ*)(mat->B)->data;
82     if (A->rmap->bs > 1) SETERRQ1(PETSC_ERR_SUP," bs=%d is not supported yet\n", A->rmap->bs);
83     nz = aa->nz + bb->nz;
84     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
85     garray = mat->garray;
86     av=aa->a; bv=bb->a;
87   }
88 
89   if (!valOnly){
90     ierr = PetscMalloc(nz*sizeof(PetscInt) ,&row);CHKERRQ(ierr);
91     ierr = PetscMalloc(nz*sizeof(PetscInt),&col);CHKERRQ(ierr);
92     ierr = PetscMalloc(nz*sizeof(PetscScalar),&val);CHKERRQ(ierr);
93     *r = row; *c = col; *v = val;
94   } else {
95     row = *r; col = *c; val = *v;
96   }
97   *nnz = nz;
98 
99   jj = 0; irow = rstart;
100   for ( i=0; i<m; i++ ) {
101     ajj = aj + ai[i];                 /* ptr to the beginning of this row */
102     countA = ai[i+1] - ai[i];
103     countB = bi[i+1] - bi[i];
104     bjj = bj + bi[i];
105 
106     /* get jB, the starting local col index for the 2nd B-part */
107     colA_start = rstart + ajj[0]; /* the smallest col index for A */
108     j=-1;
109     do {
110       j++;
111       if (j == countB) break;
112       jcol = garray[bjj[j]];
113     } while (jcol < colA_start);
114     jB = j;
115 
116     /* B-part, smaller col index */
117     colA_start = rstart + ajj[0]; /* the smallest col index for A */
118     for (j=0; j<jB; j++){
119       jcol = garray[bjj[j]];
120       if (!valOnly){
121         row[jj] = irow + shift; col[jj] = jcol + shift;
122 
123       }
124       val[jj++] = *bv++;
125     }
126     /* A-part */
127     for (j=0; j<countA; j++){
128       if (!valOnly){
129         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
130       }
131       val[jj++] = *av++;
132     }
133     /* B-part, larger col index */
134     for (j=jB; j<countB; j++){
135       if (!valOnly){
136         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
137       }
138       val[jj++] = *bv++;
139     }
140     irow++;
141   }
142 
143   PetscFunctionReturn(0);
144 }
145 
146 #undef __FUNCT__
147 #define __FUNCT__ "MatDestroy_MUMPS"
148 PetscErrorCode MatDestroy_MUMPS(Mat A)
149 {
150   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
151   PetscErrorCode ierr;
152   PetscMPIInt    size=lu->size;
153 
154   PetscFunctionBegin;
155   if (lu->CleanUpMUMPS) {
156     /* Terminate instance, deallocate memories */
157     if (size > 1){
158       ierr = PetscFree(lu->id.sol_loc);CHKERRQ(ierr);
159       ierr = VecScatterDestroy(lu->scat_rhs);CHKERRQ(ierr);
160       ierr = VecDestroy(lu->b_seq);CHKERRQ(ierr);
161       if (lu->nSolve && lu->scat_sol){ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);}
162       if (lu->nSolve && lu->x_seq){ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);}
163       ierr = PetscFree(lu->val);CHKERRQ(ierr);
164     }
165     lu->id.job=JOB_END;
166 #if defined(PETSC_USE_COMPLEX)
167     zmumps_c(&lu->id);
168 #else
169     dmumps_c(&lu->id);
170 #endif
171     ierr = PetscFree(lu->irn);CHKERRQ(ierr);
172     ierr = PetscFree(lu->jcn);CHKERRQ(ierr);
173     ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr);
174   }
175   ierr = (lu->MatDestroy)(A);CHKERRQ(ierr);
176   PetscFunctionReturn(0);
177 }
178 
179 #undef __FUNCT__
180 #define __FUNCT__ "MatSolve_MUMPS"
181 PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x)
182 {
183   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
184   PetscScalar    *array;
185   Vec            x_seq;
186   IS             is_iden,is_petsc;
187   PetscErrorCode ierr;
188   PetscInt       i;
189 
190   PetscFunctionBegin;
191   lu->id.nrhs = 1;
192   x_seq = lu->b_seq;
193   if (lu->size > 1){
194     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
195     ierr = VecScatterBegin(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
196     ierr = VecScatterEnd(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
197     if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);}
198   } else {  /* size == 1 */
199     ierr = VecCopy(b,x);CHKERRQ(ierr);
200     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
201   }
202   if (!lu->myid) { /* define rhs on the host */
203     lu->id.nrhs = 1;
204 #if defined(PETSC_USE_COMPLEX)
205     lu->id.rhs = (mumps_double_complex*)array;
206 #else
207     lu->id.rhs = array;
208 #endif
209   }
210   if (lu->size == 1){
211     ierr = VecRestoreArray(x,&array);CHKERRQ(ierr);
212   } else if (!lu->myid){
213     ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr);
214   }
215 
216   if (lu->size > 1){
217     /* distributed solution */
218     lu->id.ICNTL(21) = 1;
219     if (!lu->nSolve){
220       /* Create x_seq=sol_loc for repeated use */
221       PetscInt    lsol_loc;
222       PetscScalar *sol_loc;
223       lsol_loc = lu->id.INFO(23); /* length of sol_loc */
224       ierr = PetscMalloc((1+lsol_loc)*(sizeof(PetscScalar)+sizeof(PetscInt)),&sol_loc);CHKERRQ(ierr);
225       lu->id.isol_loc = (PetscInt *)(sol_loc + lsol_loc);
226       lu->id.lsol_loc = lsol_loc;
227 #if defined(PETSC_USE_COMPLEX)
228       lu->id.sol_loc  = (mumps_double_complex*)sol_loc;
229 #else
230       lu->id.sol_loc  = sol_loc;
231 #endif
232       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,lsol_loc,sol_loc,&lu->x_seq);CHKERRQ(ierr);
233     }
234   }
235 
236   /* solve phase */
237   /*-------------*/
238   lu->id.job = 3;
239 #if defined(PETSC_USE_COMPLEX)
240   zmumps_c(&lu->id);
241 #else
242   dmumps_c(&lu->id);
243 #endif
244   if (lu->id.INFOG(1) < 0) {
245     SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1));
246   }
247 
248   if (lu->size > 1) { /* convert mumps distributed solution to petsc mpi x */
249     if (!lu->nSolve){ /* create scatter scat_sol */
250       ierr = ISCreateStride(PETSC_COMM_SELF,lu->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */
251       for (i=0; i<lu->id.lsol_loc; i++){
252         lu->id.isol_loc[i] -= 1; /* change Fortran style to C style */
253       }
254       ierr = ISCreateGeneral(PETSC_COMM_SELF,lu->id.lsol_loc,lu->id.isol_loc,&is_petsc);CHKERRQ(ierr);  /* to */
255       ierr = VecScatterCreate(lu->x_seq,is_iden,x,is_petsc,&lu->scat_sol);CHKERRQ(ierr);
256       ierr = ISDestroy(is_iden);CHKERRQ(ierr);
257       ierr = ISDestroy(is_petsc);CHKERRQ(ierr);
258     }
259     ierr = VecScatterBegin(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
260     ierr = VecScatterEnd(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
261   }
262   lu->nSolve++;
263   PetscFunctionReturn(0);
264 }
265 
266 #if !defined(PETSC_USE_COMPLEX)
267 /*
268   input:
269    F:        numeric factor
270   output:
271    nneg:     total number of negative pivots
272    nzero:    0
273    npos:     (global dimension of F) - nneg
274 */
275 
276 #undef __FUNCT__
277 #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
278 PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
279 {
280   Mat_MUMPS      *lu =(Mat_MUMPS*)F->spptr;
281   PetscErrorCode ierr;
282   PetscMPIInt    size;
283 
284   PetscFunctionBegin;
285   ierr = MPI_Comm_size(((PetscObject)F)->comm,&size);CHKERRQ(ierr);
286   /* 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 */
287   if (size > 1 && lu->id.ICNTL(13) != 1){
288     SETERRQ1(PETSC_ERR_ARG_WRONG,"ICNTL(13)=%d. -mat_mumps_icntl_13 must be set as 1 for correct global matrix inertia\n",lu->id.INFOG(13));
289   }
290   if (nneg){
291     if (!lu->myid){
292       *nneg = lu->id.INFOG(12);
293     }
294     ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr);
295   }
296   if (nzero) *nzero = 0;
297   if (npos)  *npos  = F->rmap->N - (*nneg);
298   PetscFunctionReturn(0);
299 }
300 #endif /* !defined(PETSC_USE_COMPLEX) */
301 
302 #undef __FUNCT__
303 #define __FUNCT__ "MatFactorNumeric_MUMPS"
304 PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
305 {
306   Mat_MUMPS      *lu =(Mat_MUMPS*)(F)->spptr;
307   PetscErrorCode ierr;
308   PetscInt       rnz,nnz,nz=0,i,M=A->rmap->N,*ai,*aj,icntl;
309   PetscTruth     valOnly,flg;
310   Mat            F_diag;
311   IS             is_iden;
312   Vec            b;
313   PetscTruth     isSeqAIJ,isSeqSBAIJ;
314 
315   PetscFunctionBegin;
316   ierr = PetscTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr);
317   ierr = PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
318   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
319     (F)->ops->solve   = MatSolve_MUMPS;
320 
321     /* Initialize a MUMPS instance */
322     ierr = MPI_Comm_rank(((PetscObject)A)->comm, &lu->myid);
323     ierr = MPI_Comm_size(((PetscObject)A)->comm,&lu->size);CHKERRQ(ierr);
324     lu->id.job = JOB_INIT;
325     ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->comm_mumps));CHKERRQ(ierr);
326     lu->id.comm_fortran = MPI_Comm_c2f(lu->comm_mumps);
327 
328     /* Set mumps options */
329     ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
330     lu->id.par=1;  /* host participates factorizaton and solve */
331     lu->id.sym=lu->sym;
332     if (lu->sym == 2){
333       ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr);
334       if (flg && icntl == 1) lu->id.sym=icntl;  /* matrix is spd */
335     }
336 #if defined(PETSC_USE_COMPLEX)
337     zmumps_c(&lu->id);
338 #else
339     dmumps_c(&lu->id);
340 #endif
341 
342     if (isSeqAIJ || isSeqSBAIJ){
343       lu->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
344     } else {
345       lu->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
346     }
347 
348     icntl=-1;
349     lu->id.ICNTL(4) = 0;  /* level of printing; overwrite mumps default ICNTL(4)=2 */
350     ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
351     if ((flg && icntl > 0) || PetscLogPrintInfo) {
352       lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */
353     } else { /* no output */
354       lu->id.ICNTL(1) = 0;  /* error message, default= 6 */
355       lu->id.ICNTL(2) = 0;  /* output stream for diagnostic printing, statistics, and warning. default=0 */
356       lu->id.ICNTL(3) = 0; /* output stream for global information, default=6 */
357     }
358     ierr = PetscOptionsInt("-mat_mumps_icntl_6","ICNTL(6): column permutation and/or scaling to get a zero-free diagonal (0 to 7)","None",lu->id.ICNTL(6),&lu->id.ICNTL(6),PETSC_NULL);CHKERRQ(ierr);
359     icntl=-1;
360     ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
361     if (flg) {
362       if (icntl== 1){
363         SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n");
364       } else {
365         lu->id.ICNTL(7) = icntl;
366       }
367     }
368     ierr = PetscOptionsInt("-mat_mumps_icntl_8","ICNTL(8): scaling strategy (-2 to 7 or 77)","None",lu->id.ICNTL(8),&lu->id.ICNTL(8),PETSC_NULL);CHKERRQ(ierr);
369     ierr = PetscOptionsInt("-mat_mumps_icntl_9","ICNTL(9): A or A^T x=b to be solved. 1: A; otherwise: A^T","None",lu->id.ICNTL(9),&lu->id.ICNTL(9),PETSC_NULL);CHKERRQ(ierr);
370     ierr = PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",lu->id.ICNTL(10),&lu->id.ICNTL(10),PETSC_NULL);CHKERRQ(ierr);
371     ierr = PetscOptionsInt("-mat_mumps_icntl_11","ICNTL(11): statistics related to the linear system solved (via -ksp_view)","None",lu->id.ICNTL(11),&lu->id.ICNTL(11),PETSC_NULL);CHKERRQ(ierr);
372     ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control: defines the ordering strategy with scaling constraints (0 to 3","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr);
373     ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control: with or without ScaLAPACK","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr);
374     ierr = PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage of estimated workspace increase","None",lu->id.ICNTL(14),&lu->id.ICNTL(14),PETSC_NULL);CHKERRQ(ierr);
375     ierr = PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): Schur complement","None",lu->id.ICNTL(19),&lu->id.ICNTL(19),PETSC_NULL);CHKERRQ(ierr);
376 
377     ierr = PetscOptionsInt("-mat_mumps_icntl_22","ICNTL(22): in-core/out-of-core facility (0 or 1)","None",lu->id.ICNTL(22),&lu->id.ICNTL(22),PETSC_NULL);CHKERRQ(ierr);
378     ierr = PetscOptionsInt("-mat_mumps_icntl_23","ICNTL(23): max size of the working memory (MB) that can allocate per processor","None",lu->id.ICNTL(23),&lu->id.ICNTL(23),PETSC_NULL);CHKERRQ(ierr);
379     ierr = PetscOptionsInt("-mat_mumps_icntl_24","ICNTL(24): detection of null pivot rows (0 or 1)","None",lu->id.ICNTL(24),&lu->id.ICNTL(24),PETSC_NULL);CHKERRQ(ierr);
380     ierr = PetscOptionsInt("-mat_mumps_icntl_25","ICNTL(25): computation of a null space basis","None",lu->id.ICNTL(25),&lu->id.ICNTL(25),PETSC_NULL);CHKERRQ(ierr);
381     ierr = PetscOptionsInt("-mat_mumps_icntl_26","ICNTL(26): Schur options for right-hand side or solution vector","None",lu->id.ICNTL(26),&lu->id.ICNTL(26),PETSC_NULL);CHKERRQ(ierr);
382     ierr = PetscOptionsInt("-mat_mumps_icntl_27","ICNTL(27): experimental parameter","None",lu->id.ICNTL(27),&lu->id.ICNTL(27),PETSC_NULL);CHKERRQ(ierr);
383 
384     ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",lu->id.CNTL(1),&lu->id.CNTL(1),PETSC_NULL);CHKERRQ(ierr);
385     ierr = PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",lu->id.CNTL(2),&lu->id.CNTL(2),PETSC_NULL);CHKERRQ(ierr);
386     ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",lu->id.CNTL(3),&lu->id.CNTL(3),PETSC_NULL);CHKERRQ(ierr);
387     ierr = PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",lu->id.CNTL(4),&lu->id.CNTL(4),PETSC_NULL);CHKERRQ(ierr);
388     ierr = PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",lu->id.CNTL(5),&lu->id.CNTL(5),PETSC_NULL);CHKERRQ(ierr);
389     PetscOptionsEnd();
390   }
391 
392   /* define matrix A */
393   switch (lu->id.ICNTL(18)){
394   case 0:  /* centralized assembled matrix input (size=1) */
395     if (!lu->myid) {
396       if (isSeqAIJ){
397         Mat_SeqAIJ   *aa = (Mat_SeqAIJ*)A->data;
398         nz               = aa->nz;
399         ai = aa->i; aj = aa->j; lu->val = aa->a;
400       } else if (isSeqSBAIJ) {
401         Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data;
402         nz                  =  aa->nz;
403         ai = aa->i; aj = aa->j; lu->val = aa->a;
404       } else {
405         SETERRQ(PETSC_ERR_SUP,"No mumps factorization for this matrix type");
406       }
407       if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */
408         ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->irn);CHKERRQ(ierr);
409         ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->jcn);CHKERRQ(ierr);
410         nz = 0;
411         for (i=0; i<M; i++){
412           rnz = ai[i+1] - ai[i];
413           while (rnz--) {  /* Fortran row/col index! */
414             lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++;
415           }
416         }
417       }
418     }
419     break;
420   case 3:  /* distributed assembled matrix input (size>1) */
421     if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
422       valOnly = PETSC_FALSE;
423     } else {
424       valOnly = PETSC_TRUE; /* only update mat values, not row and col index */
425     }
426     ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
427     break;
428   default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS.");
429   }
430 
431   /* analysis phase */
432   /*----------------*/
433   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
434     lu->id.job = 1;
435 
436     lu->id.n = M;
437     switch (lu->id.ICNTL(18)){
438     case 0:  /* centralized assembled matrix input */
439       if (!lu->myid) {
440         lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
441         if (lu->id.ICNTL(6)>1){
442 #if defined(PETSC_USE_COMPLEX)
443           lu->id.a = (mumps_double_complex*)lu->val;
444 #else
445           lu->id.a = lu->val;
446 #endif
447         }
448       }
449       break;
450     case 3:  /* distributed assembled matrix input (size>1) */
451       lu->id.nz_loc = nnz;
452       lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
453       if (lu->id.ICNTL(6)>1) {
454 #if defined(PETSC_USE_COMPLEX)
455         lu->id.a_loc = (mumps_double_complex*)lu->val;
456 #else
457         lu->id.a_loc = lu->val;
458 #endif
459       }
460       /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
461       if (!lu->myid){
462         ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr);
463         ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
464       } else {
465         ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr);
466         ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
467       }
468       ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr);
469       ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
470       ierr = VecSetFromOptions(b);CHKERRQ(ierr);
471 
472       ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr);
473       ierr = ISDestroy(is_iden);CHKERRQ(ierr);
474       ierr = VecDestroy(b);CHKERRQ(ierr);
475       break;
476     }
477 #if defined(PETSC_USE_COMPLEX)
478     zmumps_c(&lu->id);
479 #else
480     dmumps_c(&lu->id);
481 #endif
482     if (lu->id.INFOG(1) < 0) {
483       SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
484     }
485   }
486 
487   /* numerical factorization phase */
488   /*-------------------------------*/
489   lu->id.job = 2;
490   if(!lu->id.ICNTL(18)) {
491     if (!lu->myid) {
492 #if defined(PETSC_USE_COMPLEX)
493       lu->id.a = (mumps_double_complex*)lu->val;
494 #else
495       lu->id.a = lu->val;
496 #endif
497     }
498   } else {
499 #if defined(PETSC_USE_COMPLEX)
500     lu->id.a_loc = (mumps_double_complex*)lu->val;
501 #else
502     lu->id.a_loc = lu->val;
503 #endif
504   }
505 #if defined(PETSC_USE_COMPLEX)
506   zmumps_c(&lu->id);
507 #else
508   dmumps_c(&lu->id);
509 #endif
510   if (lu->id.INFOG(1) < 0) {
511     if (lu->id.INFO(1) == -13) {
512       SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",lu->id.INFO(2));
513     } else {
514       SETERRQ2(PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: INFO(1)=%d, INFO(2)=%d\n",lu->id.INFO(1),lu->id.INFO(2));
515     }
516   }
517 
518   if (!lu->myid && lu->id.ICNTL(16) > 0){
519     SETERRQ1(PETSC_ERR_LIB,"  lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16));
520   }
521 
522   if (lu->size > 1){
523     if ((F)->factor == MAT_FACTOR_LU){
524       F_diag = ((Mat_MPIAIJ *)(F)->data)->A;
525     } else {
526       F_diag = ((Mat_MPISBAIJ *)(F)->data)->A;
527     }
528     F_diag->assembled = PETSC_TRUE;
529     if (lu->nSolve){
530       ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);
531       ierr = PetscFree(lu->id.sol_loc);CHKERRQ(ierr);
532       ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);
533     }
534   }
535   (F)->assembled   = PETSC_TRUE;
536   lu->matstruc      = SAME_NONZERO_PATTERN;
537   lu->CleanUpMUMPS  = PETSC_TRUE;
538   lu->nSolve        = 0;
539   PetscFunctionReturn(0);
540 }
541 
542 
543 /* Note the Petsc r and c permutations are ignored */
544 #undef __FUNCT__
545 #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
546 PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
547 {
548   Mat_MUMPS      *lu = (Mat_MUMPS*)F->spptr;
549 
550   PetscFunctionBegin;
551   lu->sym                  = 0;
552   lu->matstruc             = DIFFERENT_NONZERO_PATTERN;
553   F->ops->lufactornumeric  = MatFactorNumeric_MUMPS;
554   PetscFunctionReturn(0);
555 }
556 
557 
558 /* Note the Petsc r permutation is ignored */
559 #undef __FUNCT__
560 #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS"
561 PetscErrorCode MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
562 {
563   Mat_MUMPS      *lu = (Mat_MUMPS*)(F)->spptr;
564 
565   PetscFunctionBegin;
566   lu->sym                          = 2;
567   lu->matstruc                     = DIFFERENT_NONZERO_PATTERN;
568   (F)->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
569 #if !defined(PETSC_USE_COMPLEX)
570   (F)->ops->getinertia            = MatGetInertia_SBAIJMUMPS;
571 #endif
572   PetscFunctionReturn(0);
573 }
574 
575 #undef __FUNCT__
576 #define __FUNCT__ "MatFactorInfo_MUMPS"
577 PetscErrorCode MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) {
578   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
579   PetscErrorCode ierr;
580 
581   PetscFunctionBegin;
582   /* check if matrix is mumps type */
583   if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0);
584 
585   ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
586   ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                  %d \n",lu->id.sym);CHKERRQ(ierr);
587   ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):           %d \n",lu->id.par);CHKERRQ(ierr);
588   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):        %d \n",lu->id.ICNTL(1));CHKERRQ(ierr);
589   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg):%d \n",lu->id.ICNTL(2));CHKERRQ(ierr);
590   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):  %d \n",lu->id.ICNTL(3));CHKERRQ(ierr);
591   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):       %d \n",lu->id.ICNTL(4));CHKERRQ(ierr);
592   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):        %d \n",lu->id.ICNTL(5));CHKERRQ(ierr);
593   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):       %d \n",lu->id.ICNTL(6));CHKERRQ(ierr);
594   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (matrix ordering):         %d \n",lu->id.ICNTL(7));CHKERRQ(ierr);
595   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scalling strategy):       %d \n",lu->id.ICNTL(8));CHKERRQ(ierr);
596   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(9) (A/A^T x=b is solved):     %d \n",lu->id.ICNTL(9));CHKERRQ(ierr);
597   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr);
598   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):         %d \n",lu->id.ICNTL(11));CHKERRQ(ierr);
599   if (!lu->myid && lu->id.ICNTL(11)>0) {
600     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(4) (inf norm of input mat):        %g\n",lu->id.RINFOG(4));CHKERRQ(ierr);
601     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(5) (inf norm of solution):         %g\n",lu->id.RINFOG(5));CHKERRQ(ierr);
602     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(6) (inf norm of residual):         %g\n",lu->id.RINFOG(6));CHKERRQ(ierr);
603     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",lu->id.RINFOG(7),lu->id.RINFOG(8));CHKERRQ(ierr);
604     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(9) (error estimate):               %g \n",lu->id.RINFOG(9));CHKERRQ(ierr);
605     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",lu->id.RINFOG(10),lu->id.RINFOG(11));CHKERRQ(ierr);
606 
607   }
608   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",lu->id.ICNTL(12));CHKERRQ(ierr);
609   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",lu->id.ICNTL(13));CHKERRQ(ierr);
610   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr);
611   /* ICNTL(15-17) not used */
612   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",lu->id.ICNTL(18));CHKERRQ(ierr);
613   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Shur complement info):                       %d \n",lu->id.ICNTL(19));CHKERRQ(ierr);
614   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",lu->id.ICNTL(20));CHKERRQ(ierr);
615   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d \n",lu->id.ICNTL(21));CHKERRQ(ierr);
616   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",lu->id.ICNTL(22));CHKERRQ(ierr);
617   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",lu->id.ICNTL(23));CHKERRQ(ierr);
618 
619   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",lu->id.ICNTL(24));CHKERRQ(ierr);
620   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",lu->id.ICNTL(25));CHKERRQ(ierr);
621   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",lu->id.ICNTL(26));CHKERRQ(ierr);
622   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",lu->id.ICNTL(27));CHKERRQ(ierr);
623 
624   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",lu->id.CNTL(1));CHKERRQ(ierr);
625   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr);
626   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",lu->id.CNTL(3));CHKERRQ(ierr);
627   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",lu->id.CNTL(4));CHKERRQ(ierr);
628   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",lu->id.CNTL(5));CHKERRQ(ierr);
629 
630   /* infomation local to each processor */
631   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);}
632   ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm,"             [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr);
633   ierr = PetscSynchronizedFlush(((PetscObject)A)->comm);
634   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);}
635   ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr);
636   ierr = PetscSynchronizedFlush(((PetscObject)A)->comm);
637   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);}
638   ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr);
639   ierr = PetscSynchronizedFlush(((PetscObject)A)->comm);
640 
641   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr);}
642   ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(15));CHKERRQ(ierr);
643   ierr = PetscSynchronizedFlush(((PetscObject)A)->comm);
644 
645   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr);}
646   ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(16));CHKERRQ(ierr);
647   ierr = PetscSynchronizedFlush(((PetscObject)A)->comm);
648 
649   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);}
650   ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(23));CHKERRQ(ierr);
651   ierr = PetscSynchronizedFlush(((PetscObject)A)->comm);
652 
653   if (!lu->myid){ /* information from the host */
654     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr);
655     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr);
656     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr);
657 
658     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr);
659     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr);
660     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr);
661     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr);
662     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr);
663     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr);
664     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",lu->id.INFOG(9));CHKERRQ(ierr);
665     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",lu->id.INFOG(10));CHKERRQ(ierr);
666     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",lu->id.INFOG(11));CHKERRQ(ierr);
667     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr);
668     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr);
669     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr);
670     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr);
671     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(16) (estimated size (in MB) of all MUMPS internal data for factorization after analysis: value on the most memory consuming processor): %d \n",lu->id.INFOG(16));CHKERRQ(ierr);
672     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(17) (estimated size of all MUMPS internal data for factorization after analysis: sum over all processors): %d \n",lu->id.INFOG(17));CHKERRQ(ierr);
673     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(18) (size of all MUMPS internal data allocated during factorization: value on the most memory consuming processor): %d \n",lu->id.INFOG(18));CHKERRQ(ierr);
674     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(19) (size of all MUMPS internal data allocated during factorization: sum over all processors): %d \n",lu->id.INFOG(19));CHKERRQ(ierr);
675      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr);
676      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(21) (size in MB of memory effectively used during factorization - value on the most memory consuming processor): %d \n",lu->id.INFOG(21));CHKERRQ(ierr);
677      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(22) (size in MB of memory effectively used during factorization - sum over all processors): %d \n",lu->id.INFOG(22));CHKERRQ(ierr);
678      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",lu->id.INFOG(23));CHKERRQ(ierr);
679      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",lu->id.INFOG(24));CHKERRQ(ierr);
680      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",lu->id.INFOG(25));CHKERRQ(ierr);
681   }
682 
683   PetscFunctionReturn(0);
684 }
685 
686 #undef __FUNCT__
687 #define __FUNCT__ "MatView_MUMPS"
688 PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
689 {
690   PetscErrorCode    ierr;
691   PetscTruth        iascii;
692   PetscViewerFormat format;
693 
694   PetscFunctionBegin;
695     ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
696   if (iascii) {
697     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
698     if (format == PETSC_VIEWER_ASCII_INFO){
699       ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr);
700     }
701   }
702   PetscFunctionReturn(0);
703 }
704 
705 #undef __FUNCT__
706 #define __FUNCT__ "MatGetInfo_MUMPS"
707 PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
708 {
709     Mat_MUMPS  *lu =(Mat_MUMPS*)A->spptr;
710 
711   PetscFunctionBegin;
712   info->block_size        = 1.0;
713   info->nz_allocated      = lu->id.INFOG(20);
714   info->nz_used           = lu->id.INFOG(20);
715   info->nz_unneeded       = 0.0;
716   info->assemblies        = 0.0;
717   info->mallocs           = 0.0;
718   info->memory            = 0.0;
719   info->fill_ratio_given  = 0;
720   info->fill_ratio_needed = 0;
721   info->factor_mallocs    = 0;
722   PetscFunctionReturn(0);
723 }
724 
725 /*MC
726   MAT_SOLVER_MUMPS -  A matrix type providing direct solvers (LU and Cholesky) for
727   distributed and sequential matrices via the external package MUMPS.
728 
729   Works with MATAIJ and MATSBAIJ matrices
730 
731   Options Database Keys:
732 + -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
733 . -mat_mumps_icntl_4 <0,...,4> - print level
734 . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
735 . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
736 . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
737 . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
738 . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
739 . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
740 . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
741 . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
742 . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
743 . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
744 . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
745 - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
746 
747   Level: beginner
748 
749 .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage
750 
751 M*/
752 
753 EXTERN_C_BEGIN
754 #undef __FUNCT__
755 #define __FUNCT__ "MatFactorGetSolverPackage_mumps"
756 PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type)
757 {
758   PetscFunctionBegin;
759   *type = MAT_SOLVER_MUMPS;
760   PetscFunctionReturn(0);
761 }
762 EXTERN_C_END
763 
764 EXTERN_C_BEGIN
765 /*
766     The seq and mpi versions of this function are the same
767 */
768 #undef __FUNCT__
769 #define __FUNCT__ "MatGetFactor_seqaij_mumps"
770 PetscErrorCode MatGetFactor_seqaij_mumps(Mat A,MatFactorType ftype,Mat *F)
771 {
772   Mat            B;
773   PetscErrorCode ierr;
774   Mat_MUMPS      *mumps;
775 
776   PetscFunctionBegin;
777   if (ftype != MAT_FACTOR_LU) {
778     SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with MUMPS Cholesky, use SBAIJ matrix");
779   }
780   /* Create the factorization matrix */
781   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
782   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
783   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
784   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
785 
786   B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
787   B->ops->view             = MatView_MUMPS;
788   B->ops->getinfo          = MatGetInfo_MUMPS;
789   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
790   B->factor                = MAT_FACTOR_LU;
791 
792   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
793   mumps->CleanUpMUMPS              = PETSC_FALSE;
794   mumps->isAIJ                     = PETSC_TRUE;
795   mumps->scat_rhs                  = PETSC_NULL;
796   mumps->scat_sol                  = PETSC_NULL;
797   mumps->nSolve                    = 0;
798   mumps->MatDestroy                = B->ops->destroy;
799   B->ops->destroy                  = MatDestroy_MUMPS;
800   B->spptr                         = (void*)mumps;
801 
802   *F = B;
803   PetscFunctionReturn(0);
804 }
805 EXTERN_C_END
806 
807 EXTERN_C_BEGIN
808 #undef __FUNCT__
809 #define __FUNCT__ "MatGetFactor_mpiaij_mumps"
810 PetscErrorCode MatGetFactor_mpiaij_mumps(Mat A,MatFactorType ftype,Mat *F)
811 {
812   Mat            B;
813   PetscErrorCode ierr;
814   Mat_MUMPS      *mumps;
815 
816   PetscFunctionBegin;
817   if (ftype != MAT_FACTOR_LU) {
818     SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with MUMPS Cholesky, use SBAIJ matrix");
819   }
820   /* Create the factorization matrix */
821   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
822   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
823   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
824   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
825   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
826 
827   B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
828   B->ops->view             = MatView_MUMPS;
829   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
830   B->factor                = MAT_FACTOR_LU;
831 
832   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
833   mumps->CleanUpMUMPS              = PETSC_FALSE;
834   mumps->isAIJ                     = PETSC_TRUE;
835   mumps->scat_rhs                  = PETSC_NULL;
836   mumps->scat_sol                  = PETSC_NULL;
837   mumps->nSolve                    = 0;
838   mumps->MatDestroy                = B->ops->destroy;
839   B->ops->destroy                  = MatDestroy_MUMPS;
840   B->spptr                         = (void*)mumps;
841 
842   *F = B;
843   PetscFunctionReturn(0);
844 }
845 EXTERN_C_END
846 
847 EXTERN_C_BEGIN
848 #undef __FUNCT__
849 #define __FUNCT__ "MatGetFactor_seqsbaij_mumps"
850 PetscErrorCode MatGetFactor_seqsbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
851 {
852   Mat            B;
853   PetscErrorCode ierr;
854   Mat_MUMPS      *mumps;
855 
856   PetscFunctionBegin;
857   if (ftype != MAT_FACTOR_CHOLESKY) {
858     SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
859   }
860   /* Create the factorization matrix */
861   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
862   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
863   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
864   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
865   ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
866 
867   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS;
868   B->ops->view                   = MatView_MUMPS;
869   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
870 
871   B->factor                      = MAT_FACTOR_CHOLESKY;
872 
873   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
874   mumps->CleanUpMUMPS              = PETSC_FALSE;
875   mumps->isAIJ                     = PETSC_TRUE;
876   mumps->scat_rhs                  = PETSC_NULL;
877   mumps->scat_sol                  = PETSC_NULL;
878   mumps->nSolve                    = 0;
879   mumps->MatDestroy                = B->ops->destroy;
880   B->ops->destroy                  = MatDestroy_MUMPS;
881   B->spptr                         = (void*)mumps;
882 
883   *F = B;
884   PetscFunctionReturn(0);
885 }
886 EXTERN_C_END
887 
888 EXTERN_C_BEGIN
889 #undef __FUNCT__
890 #define __FUNCT__ "MatGetFactor_mpisbaij_mumps"
891 PetscErrorCode MatGetFactor_mpisbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
892 {
893   Mat            B;
894   PetscErrorCode ierr;
895   Mat_MUMPS      *mumps;
896 
897   PetscFunctionBegin;
898   if (ftype != MAT_FACTOR_CHOLESKY) {
899     SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
900   }
901   /* Create the factorization matrix */
902   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
903   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
904   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
905   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
906   ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
907 
908   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS;
909   B->ops->view                   = MatView_MUMPS;
910   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
911   B->factor                      = MAT_FACTOR_CHOLESKY;
912 
913   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
914   mumps->CleanUpMUMPS              = PETSC_FALSE;
915   mumps->isAIJ                     = PETSC_TRUE;
916   mumps->scat_rhs                  = PETSC_NULL;
917   mumps->scat_sol                  = PETSC_NULL;
918   mumps->nSolve                    = 0;
919   mumps->MatDestroy                = B->ops->destroy;
920   B->ops->destroy                  = MatDestroy_MUMPS;
921   B->spptr                         = (void*)mumps;
922 
923   *F = B;
924   PetscFunctionReturn(0);
925 }
926 EXTERN_C_END
927