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