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