xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision f89e2e87fbc2bf363b9a1b6d11cc8db8663cdbc0)
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"  /*I  "petscmat.h"  I*/
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 #include "../src/mat/impls/baij/seq/baij.h"
11 #include "../src/mat/impls/baij/mpi/mpibaij.h"
12 
13 EXTERN_C_BEGIN
14 #if defined(PETSC_USE_COMPLEX)
15 #include "zmumps_c.h"
16 #else
17 #include "dmumps_c.h"
18 #endif
19 EXTERN_C_END
20 #define JOB_INIT -1
21 #define JOB_END -2
22 /* macros s.t. indices match MUMPS documentation */
23 #define ICNTL(I) icntl[(I)-1]
24 #define CNTL(I) cntl[(I)-1]
25 #define INFOG(I) infog[(I)-1]
26 #define INFO(I) info[(I)-1]
27 #define RINFOG(I) rinfog[(I)-1]
28 #define RINFO(I) rinfo[(I)-1]
29 
30 typedef struct {
31 #if defined(PETSC_USE_COMPLEX)
32   ZMUMPS_STRUC_C id;
33 #else
34   DMUMPS_STRUC_C id;
35 #endif
36   MatStructure   matstruc;
37   PetscMPIInt    myid,size;
38   PetscInt       *irn,*jcn,nz,sym,nSolve;
39   PetscScalar    *val;
40   MPI_Comm       comm_mumps;
41   VecScatter     scat_rhs, scat_sol;
42   PetscTruth     isAIJ,CleanUpMUMPS,mumpsview;
43   Vec            b_seq,x_seq;
44   PetscErrorCode (*MatDestroy)(Mat);
45   PetscErrorCode (*ConvertToTriples)(Mat, int, MatReuse, int*, int**, int**, PetscScalar**);
46 } Mat_MUMPS;
47 
48 EXTERN PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);
49 
50 
51 /* MatConvertToTriples_A_B */
52 /*convert Petsc matrix to triples: row[nz], col[nz], val[nz] */
53 /*
54   input:
55     A       - matrix in aij,baij or sbaij (bs=1) format
56     shift   - 0: C style output triple; 1: Fortran style output triple.
57     reuse   - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple
58               MAT_REUSE_MATRIX:   only the values in v array are updated
59   output:
60     nnz     - dim of r, c, and v (number of local nonzero entries of A)
61     r, c, v - row and col index, matrix values (matrix triples)
62  */
63 
64 #undef __FUNCT__
65 #define __FUNCT__ "MatConvertToTriples_seqaij_seqaij"
66 PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
67 {
68   const PetscInt   *ai,*aj,*ajj,M=A->rmap->n;
69   PetscInt         nz,rnz,i,j;
70   PetscErrorCode   ierr;
71   PetscInt         *row,*col;
72   Mat_SeqAIJ       *aa=(Mat_SeqAIJ*)A->data;
73 
74   PetscFunctionBegin;
75   *v=aa->a;
76   if (reuse == MAT_INITIAL_MATRIX){
77     nz = aa->nz; ai = aa->i; aj = aa->j;
78     *nnz = nz;
79     ierr = PetscMalloc(2*nz*sizeof(PetscInt), &row);CHKERRQ(ierr);
80     col  = row + nz;
81 
82     nz = 0;
83     for(i=0; i<M; i++) {
84       rnz = ai[i+1] - ai[i];
85       ajj = aj + ai[i];
86       for(j=0; j<rnz; j++) {
87 	row[nz] = i+shift; col[nz++] = ajj[j] + shift;
88       }
89     }
90     *r = row; *c = col;
91   }
92   PetscFunctionReturn(0);
93 }
94 
95 #undef __FUNCT__
96 #define __FUNCT__ "MatConvertToTriples_seqbaij_seqaij"
97 PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
98 {
99   Mat_SeqBAIJ        *aa=(Mat_SeqBAIJ*)A->data;
100   const PetscInt     *ai,*aj,*ajj,bs=A->rmap->bs,bs2=aa->bs2,M=A->rmap->N/bs;
101   PetscInt           nz,idx=0,rnz,i,j,k,m,ii;
102   PetscErrorCode     ierr;
103   PetscInt           *row,*col;
104 
105   PetscFunctionBegin;
106   *v = aa->a;
107   if (reuse == MAT_INITIAL_MATRIX){
108     ai = aa->i; aj = aa->j;
109     nz = bs2*aa->nz;
110     *nnz = nz;
111     ierr = PetscMalloc(2*nz*sizeof(PetscInt), &row);CHKERRQ(ierr);
112     col  = row + nz;
113 
114     for(i=0; i<M; i++) {
115       ii = 0;
116       ajj = aj + ai[i];
117       rnz = ai[i+1] - ai[i];
118       for(k=0; k<rnz; k++) {
119 	for(j=0; j<bs; j++) {
120 	  for(m=0; m<bs; m++) {
121 	    row[idx]     = i*bs + m + shift;
122 	    col[idx++]   = bs*(ajj[k]) + j + shift;
123 	  }
124 	}
125       }
126     }
127     *r = row; *c = col;
128   }
129   PetscFunctionReturn(0);
130 }
131 
132 #undef __FUNCT__
133 #define __FUNCT__ "MatConvertToTriples_seqsbaij_seqsbaij"
134 PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
135 {
136   const PetscInt   *ai, *aj,*ajj,M=A->rmap->n;
137   PetscInt         nz,rnz,i,j;
138   PetscErrorCode   ierr;
139   PetscInt         *row,*col;
140   Mat_SeqSBAIJ     *aa=(Mat_SeqSBAIJ*)A->data;
141 
142   PetscFunctionBegin;
143   if (reuse == MAT_INITIAL_MATRIX){
144     nz = aa->nz;ai=aa->i; aj=aa->j;*v=aa->a;
145     *nnz = nz;
146     ierr = PetscMalloc(2*nz*sizeof(PetscInt), &row);CHKERRQ(ierr);
147     col  = row + nz;
148 
149     nz = 0;
150     for(i=0; i<M; i++) {
151       rnz = ai[i+1] - ai[i];
152       ajj = aj + ai[i];
153       for(j=0; j<rnz; j++) {
154 	row[nz] = i+shift; col[nz++] = ajj[j] + shift;
155       }
156     }
157     *r = row; *c = col;
158   }
159   PetscFunctionReturn(0);
160 }
161 
162 #undef __FUNCT__
163 #define __FUNCT__ "MatConvertToTriples_seqaij_seqsbaij"
164 PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
165 {
166   const PetscInt     *ai,*aj,*ajj,*adiag,M=A->rmap->n;
167   PetscInt           nz,rnz,i,j;
168   const PetscScalar  *av,*v1;
169   PetscScalar        *val;
170   PetscErrorCode     ierr;
171   PetscInt           *row,*col;
172   Mat_SeqSBAIJ       *aa=(Mat_SeqSBAIJ*)A->data;
173 
174   PetscFunctionBegin;
175   ai=aa->i; aj=aa->j;av=aa->a;
176   adiag=aa->diag;
177   if (reuse == MAT_INITIAL_MATRIX){
178     nz = M + (aa->nz-M)/2;
179     *nnz = nz;
180     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
181     col  = row + nz;
182     val  = (PetscScalar*)(col + nz);
183 
184     nz = 0;
185     for(i=0; i<M; i++) {
186       rnz = ai[i+1] - adiag[i];
187       ajj  = aj + adiag[i];
188       v1   = av + adiag[i];
189       for(j=0; j<rnz; j++) {
190 	row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
191       }
192     }
193     *r = row; *c = col; *v = val;
194   } else {
195     nz = 0; val = *v;
196     for(i=0; i <M; i++) {
197       rnz = ai[i+1] - adiag[i];
198       ajj = aj + adiag[i];
199       v1  = av + adiag[i];
200       for(j=0; j<rnz; j++) {
201 	val[nz++] = v1[j];
202       }
203     }
204   }
205   PetscFunctionReturn(0);
206 }
207 
208 #undef __FUNCT__
209 #define __FUNCT__ "MatConvertToTriples_mpisbaij_mpisbaij"
210 PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
211 {
212   const PetscInt     *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
213   PetscErrorCode     ierr;
214   PetscInt           rstart,nz,i,j,jj,irow,countA,countB;
215   PetscInt           *row,*col;
216   const PetscScalar  *av, *bv,*v1,*v2;
217   PetscScalar        *val;
218   Mat_MPISBAIJ       *mat =  (Mat_MPISBAIJ*)A->data;
219   Mat_SeqSBAIJ       *aa=(Mat_SeqSBAIJ*)(mat->A)->data;
220   Mat_SeqBAIJ        *bb=(Mat_SeqBAIJ*)(mat->B)->data;
221 
222   PetscFunctionBegin;
223   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
224   garray = mat->garray;
225   av=aa->a; bv=bb->a;
226 
227   if (reuse == MAT_INITIAL_MATRIX){
228     nz = aa->nz + bb->nz;
229     *nnz = nz;
230     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
231     col  = row + nz;
232     val  = (PetscScalar*)(col + nz);
233 
234     *r = row; *c = col; *v = val;
235   } else {
236     row = *r; col = *c; val = *v;
237   }
238 
239   jj = 0; irow = rstart;
240   for ( i=0; i<m; i++ ) {
241     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
242     countA = ai[i+1] - ai[i];
243     countB = bi[i+1] - bi[i];
244     bjj    = bj + bi[i];
245     v1     = av + ai[i];
246     v2     = bv + bi[i];
247 
248     /* A-part */
249     for (j=0; j<countA; j++){
250       if (reuse == MAT_INITIAL_MATRIX) {
251         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
252       }
253       val[jj++] = v1[j];
254     }
255 
256     /* B-part */
257     for(j=0; j < countB; j++){
258       if (reuse == MAT_INITIAL_MATRIX) {
259 	row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
260       }
261       val[jj++] = v2[j];
262     }
263     irow++;
264   }
265   PetscFunctionReturn(0);
266 }
267 
268 #undef __FUNCT__
269 #define __FUNCT__ "MatConvertToTriples_mpiaij_mpiaij"
270 PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
271 {
272   const PetscInt     *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
273   PetscErrorCode     ierr;
274   PetscInt           rstart,nz,i,j,jj,irow,countA,countB;
275   PetscInt           *row,*col;
276   const PetscScalar  *av, *bv,*v1,*v2;
277   PetscScalar        *val;
278   Mat_MPIAIJ         *mat =  (Mat_MPIAIJ*)A->data;
279   Mat_SeqAIJ         *aa=(Mat_SeqAIJ*)(mat->A)->data;
280   Mat_SeqAIJ         *bb=(Mat_SeqAIJ*)(mat->B)->data;
281 
282   PetscFunctionBegin;
283   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
284   garray = mat->garray;
285   av=aa->a; bv=bb->a;
286 
287   if (reuse == MAT_INITIAL_MATRIX){
288     nz = aa->nz + bb->nz;
289     *nnz = nz;
290     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
291     col  = row + nz;
292     val  = (PetscScalar*)(col + nz);
293 
294     *r = row; *c = col; *v = val;
295   } else {
296     row = *r; col = *c; val = *v;
297   }
298 
299   jj = 0; irow = rstart;
300   for ( i=0; i<m; i++ ) {
301     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
302     countA = ai[i+1] - ai[i];
303     countB = bi[i+1] - bi[i];
304     bjj    = bj + bi[i];
305     v1     = av + ai[i];
306     v2     = bv + bi[i];
307 
308     /* A-part */
309     for (j=0; j<countA; j++){
310       if (reuse == MAT_INITIAL_MATRIX){
311         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
312       }
313       val[jj++] = v1[j];
314     }
315 
316     /* B-part */
317     for(j=0; j < countB; j++){
318       if (reuse == MAT_INITIAL_MATRIX){
319 	row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
320       }
321       val[jj++] = v2[j];
322     }
323     irow++;
324   }
325   PetscFunctionReturn(0);
326 }
327 
328 #undef __FUNCT__
329 #define __FUNCT__ "MatConvertToTriples_mpibaij_mpiaij"
330 PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
331 {
332   Mat_MPIBAIJ        *mat =  (Mat_MPIBAIJ*)A->data;
333   Mat_SeqBAIJ        *aa=(Mat_SeqBAIJ*)(mat->A)->data;
334   Mat_SeqBAIJ        *bb=(Mat_SeqBAIJ*)(mat->B)->data;
335   const PetscInt     *ai = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj;
336   const PetscInt     *garray = mat->garray,mbs=mat->mbs,rstartbs=mat->rstartbs;
337   const PetscInt     bs = A->rmap->bs,bs2=mat->bs2;
338   PetscErrorCode     ierr;
339   PetscInt           nz,i,j,k,n,jj,irow,countA,countB,idx;
340   PetscInt           *row,*col;
341   const PetscScalar  *av=aa->a, *bv=bb->a,*v1,*v2;
342   PetscScalar        *val;
343 
344   PetscFunctionBegin;
345 
346   if (reuse == MAT_INITIAL_MATRIX) {
347     nz = bs2*(aa->nz + bb->nz);
348     *nnz = nz;
349     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
350     col  = row + nz;
351     val  = (PetscScalar*)(col + nz);
352 
353     *r = row; *c = col; *v = val;
354   } else {
355     row = *r; col = *c; val = *v;
356   }
357 
358   jj = 0; irow = rstartbs;
359   for ( i=0; i<mbs; i++ ) {
360     countA = ai[i+1] - ai[i];
361     countB = bi[i+1] - bi[i];
362     ajj    = aj + ai[i];
363     bjj    = bj + bi[i];
364     v1     = av + bs2*ai[i];
365     v2     = bv + bs2*bi[i];
366 
367     idx = 0;
368     /* A-part */
369     for (k=0; k<countA; k++){
370       for (j=0; j<bs; j++) {
371 	for (n=0; n<bs; n++) {
372 	  if (reuse == MAT_INITIAL_MATRIX){
373 	    row[jj] = bs*irow + n + shift;
374 	    col[jj] = bs*(rstartbs + ajj[k]) + j + shift;
375 	  }
376 	  val[jj++] = v1[idx++];
377 	}
378       }
379     }
380 
381     idx = 0;
382     /* B-part */
383     for(k=0; k<countB; k++){
384       for (j=0; j<bs; j++) {
385 	for (n=0; n<bs; n++) {
386 	  if (reuse == MAT_INITIAL_MATRIX){
387 	    row[jj] = bs*irow + n + shift;
388 	    col[jj] = bs*(garray[bjj[k]]) + j + shift;
389 	  }
390 	  val[jj++] = bv[idx++];
391 	}
392       }
393     }
394     irow++;
395   }
396   PetscFunctionReturn(0);
397 }
398 
399 #undef __FUNCT__
400 #define __FUNCT__ "MatConvertToTriples_mpiaij_mpisbaij"
401 PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
402 {
403   const PetscInt     *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
404   PetscErrorCode     ierr;
405   PetscInt           rstart,nz,nza,nzb_low,i,j,jj,irow,countA,countB;
406   PetscInt           *row,*col;
407   const PetscScalar  *av, *bv,*v1,*v2;
408   PetscScalar        *val;
409   Mat_MPIAIJ         *mat =  (Mat_MPIAIJ*)A->data;
410   Mat_SeqAIJ         *aa=(Mat_SeqAIJ*)(mat->A)->data;
411   Mat_SeqAIJ         *bb=(Mat_SeqAIJ*)(mat->B)->data;
412 
413   PetscFunctionBegin;
414   ai=aa->i; aj=aa->j; adiag=aa->diag;
415   bi=bb->i; bj=bb->j; garray = mat->garray;
416   av=aa->a; bv=bb->a;
417   rstart = A->rmap->rstart;
418 
419   if (reuse == MAT_INITIAL_MATRIX) {
420     nza = 0;nzb_low = 0;
421     for(i=0; i<m; i++){
422       nza     = nza + (ai[i+1] - adiag[i]);
423       countB  = bi[i+1] - bi[i];
424       bjj     = bj + bi[i];
425 
426       j = 0;
427       while(garray[bjj[j]] < rstart) {
428 	if(j == countB) break;
429 	j++;nzb_low++;
430       }
431     }
432     /* Total nz = nz for the upper triangular A part + nz for the 2nd B part */
433     nz = nza + (bb->nz - nzb_low);
434     *nnz = nz;
435     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
436     col  = row + nz;
437     val  = (PetscScalar*)(col + nz);
438 
439     *r = row; *c = col; *v = val;
440   } else {
441     row = *r; col = *c; val = *v;
442   }
443 
444   jj = 0; irow = rstart;
445   for ( i=0; i<m; i++ ) {
446     ajj    = aj + adiag[i];                 /* ptr to the beginning of the diagonal of this row */
447     v1     = av + adiag[i];
448     countA = ai[i+1] - adiag[i];
449     countB = bi[i+1] - bi[i];
450     bjj    = bj + bi[i];
451     v2     = bv + bi[i];
452 
453      /* A-part */
454     for (j=0; j<countA; j++){
455       if (reuse == MAT_INITIAL_MATRIX) {
456         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
457       }
458       val[jj++] = v1[j];
459     }
460 
461     /* B-part */
462     for(j=0; j < countB; j++){
463       if (garray[bjj[j]] > rstart) {
464 	if (reuse == MAT_INITIAL_MATRIX) {
465 	  row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
466 	}
467 	val[jj++] = v2[j];
468       }
469     }
470     irow++;
471   }
472   PetscFunctionReturn(0);
473 }
474 
475 #undef __FUNCT__
476 #define __FUNCT__ "MatDestroy_MUMPS"
477 PetscErrorCode MatDestroy_MUMPS(Mat A)
478 {
479   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
480   PetscErrorCode ierr;
481   PetscMPIInt    size=lu->size;
482 
483   PetscFunctionBegin;
484   if (lu->CleanUpMUMPS) {
485     /* Terminate instance, deallocate memories */
486     if (size > 1){
487       ierr = PetscFree2(lu->id.sol_loc,lu->id.isol_loc);CHKERRQ(ierr);
488       ierr = VecScatterDestroy(lu->scat_rhs);CHKERRQ(ierr);
489       ierr = VecDestroy(lu->b_seq);CHKERRQ(ierr);
490       if (lu->nSolve && lu->scat_sol){ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);}
491       if (lu->nSolve && lu->x_seq){ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);}
492     }
493     ierr = PetscFree(lu->irn);CHKERRQ(ierr);
494     lu->id.job=JOB_END;
495 #if defined(PETSC_USE_COMPLEX)
496     zmumps_c(&lu->id);
497 #else
498     dmumps_c(&lu->id);
499 #endif
500     ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr);
501   }
502   /* clear composed functions */
503   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatFactorGetSolverPackage_C","",PETSC_NULL);CHKERRQ(ierr);
504   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMumpsSetIcntl_C","",PETSC_NULL);CHKERRQ(ierr);
505   ierr = (lu->MatDestroy)(A);CHKERRQ(ierr);
506   PetscFunctionReturn(0);
507 }
508 
509 #undef __FUNCT__
510 #define __FUNCT__ "MatSolve_MUMPS"
511 PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x)
512 {
513   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
514   PetscScalar    *array;
515   Vec            b_seq;
516   IS             is_iden,is_petsc;
517   PetscErrorCode ierr;
518   PetscInt       i;
519 
520   PetscFunctionBegin;
521   lu->id.nrhs = 1;
522   b_seq = lu->b_seq;
523   if (lu->size > 1){
524     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
525     ierr = VecScatterBegin(lu->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
526     ierr = VecScatterEnd(lu->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
527     if (!lu->myid) {ierr = VecGetArray(b_seq,&array);CHKERRQ(ierr);}
528   } else {  /* size == 1 */
529     ierr = VecCopy(b,x);CHKERRQ(ierr);
530     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
531   }
532   if (!lu->myid) { /* define rhs on the host */
533     lu->id.nrhs = 1;
534 #if defined(PETSC_USE_COMPLEX)
535     lu->id.rhs = (mumps_double_complex*)array;
536 #else
537     lu->id.rhs = array;
538 #endif
539   }
540 
541   /* solve phase */
542   /*-------------*/
543   lu->id.job = 3;
544 #if defined(PETSC_USE_COMPLEX)
545   zmumps_c(&lu->id);
546 #else
547   dmumps_c(&lu->id);
548 #endif
549   if (lu->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1));
550 
551   if (lu->size > 1) { /* convert mumps distributed solution to petsc mpi x */
552     if (!lu->nSolve){ /* create scatter scat_sol */
553       ierr = ISCreateStride(PETSC_COMM_SELF,lu->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */
554       for (i=0; i<lu->id.lsol_loc; i++){
555         lu->id.isol_loc[i] -= 1; /* change Fortran style to C style */
556       }
557       ierr = ISCreateGeneral(PETSC_COMM_SELF,lu->id.lsol_loc,lu->id.isol_loc,&is_petsc);CHKERRQ(ierr);  /* to */
558       ierr = VecScatterCreate(lu->x_seq,is_iden,x,is_petsc,&lu->scat_sol);CHKERRQ(ierr);
559       ierr = ISDestroy(is_iden);CHKERRQ(ierr);
560       ierr = ISDestroy(is_petsc);CHKERRQ(ierr);
561     }
562     ierr = VecScatterBegin(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
563     ierr = VecScatterEnd(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
564   }
565   lu->nSolve++;
566   PetscFunctionReturn(0);
567 }
568 
569 #if !defined(PETSC_USE_COMPLEX)
570 /*
571   input:
572    F:        numeric factor
573   output:
574    nneg:     total number of negative pivots
575    nzero:    0
576    npos:     (global dimension of F) - nneg
577 */
578 
579 #undef __FUNCT__
580 #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
581 PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
582 {
583   Mat_MUMPS      *lu =(Mat_MUMPS*)F->spptr;
584   PetscErrorCode ierr;
585   PetscMPIInt    size;
586 
587   PetscFunctionBegin;
588   ierr = MPI_Comm_size(((PetscObject)F)->comm,&size);CHKERRQ(ierr);
589   /* 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 */
590   if (size > 1 && lu->id.ICNTL(13) != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"ICNTL(13)=%d. -mat_mumps_icntl_13 must be set as 1 for correct global matrix inertia\n",lu->id.INFOG(13));
591   if (nneg){
592     if (!lu->myid){
593       *nneg = lu->id.INFOG(12);
594     }
595     ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr);
596   }
597   if (nzero) *nzero = 0;
598   if (npos)  *npos  = F->rmap->N - (*nneg);
599   PetscFunctionReturn(0);
600 }
601 #endif /* !defined(PETSC_USE_COMPLEX) */
602 
603 #undef __FUNCT__
604 #define __FUNCT__ "MatFactorNumeric_MUMPS"
605 PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
606 {
607   Mat_MUMPS       *lu =(Mat_MUMPS*)(F)->spptr;
608   PetscErrorCode  ierr;
609   MatReuse        reuse;
610   Mat             F_diag;
611   PetscTruth      isMPIAIJ;
612 
613   PetscFunctionBegin;
614   reuse = MAT_REUSE_MATRIX;
615   ierr = (*lu->ConvertToTriples)(A, 1, reuse, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
616 
617   /* numerical factorization phase */
618   /*-------------------------------*/
619   lu->id.job = 2;
620   if(!lu->id.ICNTL(18)) {
621     if (!lu->myid) {
622 #if defined(PETSC_USE_COMPLEX)
623       lu->id.a = (mumps_double_complex*)lu->val;
624 #else
625       lu->id.a = lu->val;
626 #endif
627     }
628   } else {
629 #if defined(PETSC_USE_COMPLEX)
630     lu->id.a_loc = (mumps_double_complex*)lu->val;
631 #else
632     lu->id.a_loc = lu->val;
633 #endif
634   }
635 #if defined(PETSC_USE_COMPLEX)
636   zmumps_c(&lu->id);
637 #else
638   dmumps_c(&lu->id);
639 #endif
640   if (lu->id.INFOG(1) < 0) {
641     if (lu->id.INFO(1) == -13) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",lu->id.INFO(2));
642     else SETERRQ2(PETSC_COMM_SELF,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));
643   }
644   if (!lu->myid && lu->id.ICNTL(16) > 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"  lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16));
645 
646   if (lu->size > 1){
647     ierr = PetscTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr);
648     if(isMPIAIJ) {
649       F_diag = ((Mat_MPIAIJ *)(F)->data)->A;
650     } else {
651       F_diag = ((Mat_MPISBAIJ *)(F)->data)->A;
652     }
653     F_diag->assembled = PETSC_TRUE;
654     if (lu->nSolve){
655       ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);
656       ierr = PetscFree2(lu->id.sol_loc,lu->id.isol_loc);CHKERRQ(ierr);
657       ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);
658     }
659   }
660   (F)->assembled   = PETSC_TRUE;
661   lu->matstruc     = SAME_NONZERO_PATTERN;
662   lu->CleanUpMUMPS = PETSC_TRUE;
663   lu->nSolve       = 0;
664 
665   if (lu->size > 1){
666     /* distributed solution */
667     lu->id.ICNTL(21) = 1;
668     if (!lu->nSolve){
669       /* Create x_seq=sol_loc for repeated use */
670       PetscInt    lsol_loc;
671       PetscScalar *sol_loc;
672       lsol_loc = lu->id.INFO(23); /* length of sol_loc */
673       ierr = PetscMalloc2(lsol_loc,PetscScalar,&sol_loc,lsol_loc,PetscInt,&lu->id.isol_loc);CHKERRQ(ierr);
674       lu->id.lsol_loc = lsol_loc;
675 #if defined(PETSC_USE_COMPLEX)
676       lu->id.sol_loc  = (mumps_double_complex*)sol_loc;
677 #else
678       lu->id.sol_loc  = sol_loc;
679 #endif
680       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,lsol_loc,sol_loc,&lu->x_seq);CHKERRQ(ierr);
681     }
682   }
683   PetscFunctionReturn(0);
684 }
685 
686 #undef __FUNCT__
687 #define __FUNCT__ "PetscSetMUMPSOptions"
688 PetscErrorCode PetscSetMUMPSOptions(Mat F, Mat A)
689 {
690   Mat_MUMPS        *lu = (Mat_MUMPS*)F->spptr;
691   PetscErrorCode   ierr;
692   PetscInt         icntl;
693   PetscTruth       flg;
694 
695   PetscFunctionBegin;
696   ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
697   ierr = PetscOptionsTruth("-mat_mumps_view","View MUMPS parameters","None",lu->mumpsview,&lu->mumpsview,PETSC_NULL);CHKERRQ(ierr);
698   if (lu->size == 1){
699     lu->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
700   } else {
701     lu->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
702   }
703 
704   icntl=-1;
705   lu->id.ICNTL(4) = 0;  /* level of printing; overwrite mumps default ICNTL(4)=2 */
706   ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
707   if ((flg && icntl > 0) || PetscLogPrintInfo) {
708     lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */
709   } else { /* no output */
710     lu->id.ICNTL(1) = 0;  /* error message, default= 6 */
711     lu->id.ICNTL(2) = 0;  /* output stream for diagnostic printing, statistics, and warning. default=0 */
712     lu->id.ICNTL(3) = 0; /* output stream for global information, default=6 */
713   }
714   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);
715   icntl=-1;
716   ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
717   if (flg) {
718     if (icntl== 1){
719       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n");
720     } else {
721       lu->id.ICNTL(7) = icntl;
722     }
723   }
724   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);
725   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);
726   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);
727   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);
728   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);
729   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);
730   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);
731   ierr = PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): Schur complement","None",lu->id.ICNTL(19),&lu->id.ICNTL(19),PETSC_NULL);CHKERRQ(ierr);
732 
733   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);
734   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);
735   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);
736   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);
737   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);
738   ierr = PetscOptionsInt("-mat_mumps_icntl_27","ICNTL(27): experimental parameter","None",lu->id.ICNTL(27),&lu->id.ICNTL(27),PETSC_NULL);CHKERRQ(ierr);
739 
740   ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",lu->id.CNTL(1),&lu->id.CNTL(1),PETSC_NULL);CHKERRQ(ierr);
741   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);
742   ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",lu->id.CNTL(3),&lu->id.CNTL(3),PETSC_NULL);CHKERRQ(ierr);
743   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);
744   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);
745   PetscOptionsEnd();
746   PetscFunctionReturn(0);
747 }
748 
749 #undef __FUNCT__
750 #define __FUNCT__ "PetscInitializeMUMPS"
751 PetscErrorCode PetscInitializeMUMPS(Mat F)
752 {
753   Mat_MUMPS       *lu = (Mat_MUMPS*)F->spptr;
754   PetscErrorCode  ierr;
755   PetscInt        icntl;
756   PetscTruth      flg;
757 
758   PetscFunctionBegin;
759   lu->id.job = JOB_INIT;
760   lu->id.par=1;  /* host participates factorizaton and solve */
761   lu->id.sym=lu->sym;
762   if (lu->sym == 2){
763     ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr);
764     if (flg && icntl == 1) lu->id.sym=icntl;  /* matrix is spd */
765   }
766 #if defined(PETSC_USE_COMPLEX)
767   zmumps_c(&lu->id);
768 #else
769   dmumps_c(&lu->id);
770 #endif
771   PetscFunctionReturn(0);
772 }
773 
774 /* Note the Petsc r and c permutations are ignored */
775 #undef __FUNCT__
776 #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
777 PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
778 {
779   Mat_MUMPS          *lu = (Mat_MUMPS*)F->spptr;
780   PetscErrorCode     ierr;
781   MatReuse           reuse;
782   Vec                b;
783   IS                 is_iden;
784   const PetscInt     M = A->rmap->N;
785 
786   PetscFunctionBegin;
787   lu->sym                  = 0;
788   lu->matstruc             = DIFFERENT_NONZERO_PATTERN;
789 
790   ierr = MPI_Comm_rank(((PetscObject)A)->comm, &lu->myid);
791   ierr = MPI_Comm_size(((PetscObject)A)->comm,&lu->size);CHKERRQ(ierr);
792   ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->comm_mumps));CHKERRQ(ierr);
793   lu->id.comm_fortran = MPI_Comm_c2f(lu->comm_mumps);
794 
795   /* Initialize a MUMPS instance */
796   ierr = PetscInitializeMUMPS(F);CHKERRQ(ierr);
797   /* Set MUMPS options */
798   ierr = PetscSetMUMPSOptions(F,A);CHKERRQ(ierr);
799 
800   reuse = MAT_INITIAL_MATRIX;
801   ierr = (*lu->ConvertToTriples)(A, 1, reuse, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
802 
803   /* analysis phase */
804   /*----------------*/
805   lu->id.job = 1;
806   lu->id.n = M;
807   switch (lu->id.ICNTL(18)){
808   case 0:  /* centralized assembled matrix input */
809     if (!lu->myid) {
810       lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
811       if (lu->id.ICNTL(6)>1){
812 #if defined(PETSC_USE_COMPLEX)
813         lu->id.a = (mumps_double_complex*)lu->val;
814 #else
815         lu->id.a = lu->val;
816 #endif
817       }
818     }
819     break;
820   case 3:  /* distributed assembled matrix input (size>1) */
821     lu->id.nz_loc = lu->nz;
822     lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
823     if (lu->id.ICNTL(6)>1) {
824 #if defined(PETSC_USE_COMPLEX)
825       lu->id.a_loc = (mumps_double_complex*)lu->val;
826 #else
827       lu->id.a_loc = lu->val;
828 #endif
829     }
830     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
831     if (!lu->myid){
832       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr);
833       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
834     } else {
835       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr);
836       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
837     }
838     ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr);
839     ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
840     ierr = VecSetFromOptions(b);CHKERRQ(ierr);
841 
842     ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr);
843     ierr = ISDestroy(is_iden);CHKERRQ(ierr);
844     ierr = VecDestroy(b);CHKERRQ(ierr);
845     break;
846     }
847 #if defined(PETSC_USE_COMPLEX)
848   zmumps_c(&lu->id);
849 #else
850   dmumps_c(&lu->id);
851 #endif
852   if (lu->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
853 
854   F->ops->lufactornumeric  = MatFactorNumeric_MUMPS;
855   F->ops->solve            = MatSolve_MUMPS;
856   PetscFunctionReturn(0);
857 }
858 
859 /* Note the Petsc r and c permutations are ignored */
860 #undef __FUNCT__
861 #define __FUNCT__ "MatLUFactorSymbolic_BAIJMUMPS"
862 PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
863 {
864 
865   Mat_MUMPS       *lu = (Mat_MUMPS*)F->spptr;
866   PetscErrorCode  ierr;
867   MatReuse        reuse;
868   Vec             b;
869   IS              is_iden;
870   const PetscInt  M = A->rmap->N;
871 
872   PetscFunctionBegin;
873   lu->sym                  = 0;
874   lu->matstruc             = DIFFERENT_NONZERO_PATTERN;
875   ierr = MPI_Comm_rank(((PetscObject)A)->comm, &lu->myid);
876   ierr = MPI_Comm_size(((PetscObject)A)->comm,&lu->size);CHKERRQ(ierr);
877   ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->comm_mumps));CHKERRQ(ierr);
878   lu->id.comm_fortran = MPI_Comm_c2f(lu->comm_mumps);
879 
880   /* Initialize a MUMPS instance */
881   ierr = PetscInitializeMUMPS(F);CHKERRQ(ierr);
882   /* Set MUMPS options */
883   ierr = PetscSetMUMPSOptions(F,A);CHKERRQ(ierr);
884 
885   reuse = MAT_INITIAL_MATRIX;
886   ierr = (*lu->ConvertToTriples)(A, 1, reuse, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
887 
888   /* analysis phase */
889   /*----------------*/
890   lu->id.job = 1;
891   lu->id.n = M;
892   switch (lu->id.ICNTL(18)){
893   case 0:  /* centralized assembled matrix input */
894     if (!lu->myid) {
895       lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
896       if (lu->id.ICNTL(6)>1){
897 #if defined(PETSC_USE_COMPLEX)
898         lu->id.a = (mumps_double_complex*)lu->val;
899 #else
900         lu->id.a = lu->val;
901 #endif
902       }
903     }
904     break;
905   case 3:  /* distributed assembled matrix input (size>1) */
906     lu->id.nz_loc = lu->nz;
907     lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
908     if (lu->id.ICNTL(6)>1) {
909 #if defined(PETSC_USE_COMPLEX)
910       lu->id.a_loc = (mumps_double_complex*)lu->val;
911 #else
912       lu->id.a_loc = lu->val;
913 #endif
914     }
915     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
916     if (!lu->myid){
917       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr);
918       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
919     } else {
920       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr);
921       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
922     }
923     ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr);
924     ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
925     ierr = VecSetFromOptions(b);CHKERRQ(ierr);
926 
927     ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr);
928     ierr = ISDestroy(is_iden);CHKERRQ(ierr);
929     ierr = VecDestroy(b);CHKERRQ(ierr);
930     break;
931     }
932 #if defined(PETSC_USE_COMPLEX)
933   zmumps_c(&lu->id);
934 #else
935   dmumps_c(&lu->id);
936 #endif
937   if (lu->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
938 
939   F->ops->lufactornumeric  = MatFactorNumeric_MUMPS;
940   F->ops->solve            = MatSolve_MUMPS;
941   PetscFunctionReturn(0);
942 }
943 
944 /* Note the Petsc r permutation is ignored */
945 #undef __FUNCT__
946 #define __FUNCT__ "MatCholeskyFactorSymbolic_MUMPS"
947 PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
948 {
949   Mat_MUMPS          *lu = (Mat_MUMPS*)F->spptr;
950   PetscErrorCode     ierr;
951   MatReuse           reuse;
952   Vec                b;
953   IS                 is_iden;
954   const PetscInt     M = A->rmap->N;
955 
956   PetscFunctionBegin;
957   lu->sym                          = 2;
958   lu->matstruc                     = DIFFERENT_NONZERO_PATTERN;
959   ierr = MPI_Comm_rank(((PetscObject)A)->comm, &lu->myid);
960   ierr = MPI_Comm_size(((PetscObject)A)->comm,&lu->size);CHKERRQ(ierr);
961   ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->comm_mumps));CHKERRQ(ierr);
962   lu->id.comm_fortran = MPI_Comm_c2f(lu->comm_mumps);
963 
964   /* Initialize a MUMPS instance */
965   ierr = PetscInitializeMUMPS(F);CHKERRQ(ierr);
966   /* Set MUMPS options */
967   ierr = PetscSetMUMPSOptions(F,A);CHKERRQ(ierr);
968 
969   reuse = MAT_INITIAL_MATRIX;
970   ierr = (*lu->ConvertToTriples)(A, 1 , reuse, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
971 
972   /* analysis phase */
973   /*----------------*/
974   lu->id.job = 1;
975   lu->id.n = M;
976   switch (lu->id.ICNTL(18)){
977   case 0:  /* centralized assembled matrix input */
978     if (!lu->myid) {
979       lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
980       if (lu->id.ICNTL(6)>1){
981 #if defined(PETSC_USE_COMPLEX)
982         lu->id.a = (mumps_double_complex*)lu->val;
983 #else
984         lu->id.a = lu->val;
985 #endif
986       }
987     }
988     break;
989   case 3:  /* distributed assembled matrix input (size>1) */
990     lu->id.nz_loc = lu->nz;
991     lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
992     if (lu->id.ICNTL(6)>1) {
993 #if defined(PETSC_USE_COMPLEX)
994       lu->id.a_loc = (mumps_double_complex*)lu->val;
995 #else
996       lu->id.a_loc = lu->val;
997 #endif
998     }
999     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1000     if (!lu->myid){
1001       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr);
1002       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
1003     } else {
1004       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr);
1005       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
1006     }
1007     ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr);
1008     ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
1009     ierr = VecSetFromOptions(b);CHKERRQ(ierr);
1010 
1011     ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr);
1012     ierr = ISDestroy(is_iden);CHKERRQ(ierr);
1013     ierr = VecDestroy(b);CHKERRQ(ierr);
1014     break;
1015     }
1016 #if defined(PETSC_USE_COMPLEX)
1017   zmumps_c(&lu->id);
1018 #else
1019   dmumps_c(&lu->id);
1020 #endif
1021   if (lu->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
1022 
1023   F->ops->choleskyfactornumeric =  MatFactorNumeric_MUMPS;
1024   F->ops->solve                 =  MatSolve_MUMPS;
1025 #if !defined(PETSC_USE_COMPLEX)
1026   (F)->ops->getinertia          =  MatGetInertia_SBAIJMUMPS;
1027 #endif
1028   PetscFunctionReturn(0);
1029 }
1030 
1031 #undef __FUNCT__
1032 #define __FUNCT__ "MatFactorInfo_MUMPS"
1033 PetscErrorCode MatFactorInfo_MUMPS(Mat A,PetscViewer viewer)
1034 {
1035   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
1036   PetscErrorCode ierr;
1037 
1038   PetscFunctionBegin;
1039   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):        %d \n",lu->id.ICNTL(1));CHKERRQ(ierr);
1040   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg):%d \n",lu->id.ICNTL(2));CHKERRQ(ierr);
1041   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):  %d \n",lu->id.ICNTL(3));CHKERRQ(ierr);
1042   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):       %d \n",lu->id.ICNTL(4));CHKERRQ(ierr);
1043   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):        %d \n",lu->id.ICNTL(5));CHKERRQ(ierr);
1044   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):       %d \n",lu->id.ICNTL(6));CHKERRQ(ierr);
1045   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (matrix ordering):         %d \n",lu->id.ICNTL(7));CHKERRQ(ierr);
1046   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scalling strategy):       %d \n",lu->id.ICNTL(8));CHKERRQ(ierr);
1047   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(9) (A/A^T x=b is solved):     %d \n",lu->id.ICNTL(9));CHKERRQ(ierr);
1048   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr);
1049   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):         %d \n",lu->id.ICNTL(11));CHKERRQ(ierr);
1050   if (lu->id.ICNTL(11)>0) {
1051     ierr = PetscViewerASCIIPrintf(viewer,"        RINFOG(4) (inf norm of input mat):        %g\n",lu->id.RINFOG(4));CHKERRQ(ierr);
1052     ierr = PetscViewerASCIIPrintf(viewer,"        RINFOG(5) (inf norm of solution):         %g\n",lu->id.RINFOG(5));CHKERRQ(ierr);
1053     ierr = PetscViewerASCIIPrintf(viewer,"        RINFOG(6) (inf norm of residual):         %g\n",lu->id.RINFOG(6));CHKERRQ(ierr);
1054     ierr = PetscViewerASCIIPrintf(viewer,"        RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",lu->id.RINFOG(7),lu->id.RINFOG(8));CHKERRQ(ierr);
1055     ierr = PetscViewerASCIIPrintf(viewer,"        RINFOG(9) (error estimate):               %g \n",lu->id.RINFOG(9));CHKERRQ(ierr);
1056     ierr = PetscViewerASCIIPrintf(viewer,"        RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",lu->id.RINFOG(10),lu->id.RINFOG(11));CHKERRQ(ierr);
1057 
1058   }
1059   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",lu->id.ICNTL(12));CHKERRQ(ierr);
1060   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",lu->id.ICNTL(13));CHKERRQ(ierr);
1061   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr);
1062   /* ICNTL(15-17) not used */
1063   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",lu->id.ICNTL(18));CHKERRQ(ierr);
1064   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Shur complement info):                       %d \n",lu->id.ICNTL(19));CHKERRQ(ierr);
1065   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",lu->id.ICNTL(20));CHKERRQ(ierr);
1066   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d \n",lu->id.ICNTL(21));CHKERRQ(ierr);
1067   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",lu->id.ICNTL(22));CHKERRQ(ierr);
1068   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",lu->id.ICNTL(23));CHKERRQ(ierr);
1069 
1070   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",lu->id.ICNTL(24));CHKERRQ(ierr);
1071   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",lu->id.ICNTL(25));CHKERRQ(ierr);
1072   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",lu->id.ICNTL(26));CHKERRQ(ierr);
1073   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",lu->id.ICNTL(27));CHKERRQ(ierr);
1074 
1075   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",lu->id.CNTL(1));CHKERRQ(ierr);
1076   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr);
1077   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",lu->id.CNTL(3));CHKERRQ(ierr);
1078   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",lu->id.CNTL(4));CHKERRQ(ierr);
1079   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",lu->id.CNTL(5));CHKERRQ(ierr);
1080 
1081   /* infomation local to each processor */
1082   ierr = PetscViewerASCIIPrintf(viewer, "              RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);
1083   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"              [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr);
1084   ierr = PetscViewerFlush(viewer);
1085   ierr = PetscViewerASCIIPrintf(viewer, "              RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);
1086   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"              [%d]  %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr);
1087   ierr = PetscViewerFlush(viewer);
1088   ierr = PetscViewerASCIIPrintf(viewer, "              RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);
1089   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"              [%d]  %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr);
1090   ierr = PetscViewerFlush(viewer);
1091 
1092   ierr = PetscViewerASCIIPrintf(viewer, "              INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr);
1093   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"              [%d] %d \n",lu->myid,lu->id.INFO(15));CHKERRQ(ierr);
1094   ierr = PetscViewerFlush(viewer);
1095 
1096   ierr = PetscViewerASCIIPrintf(viewer, "              INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr);
1097   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"              [%d] %d \n",lu->myid,lu->id.INFO(16));CHKERRQ(ierr);
1098   ierr = PetscViewerFlush(viewer);
1099 
1100   ierr = PetscViewerASCIIPrintf(viewer, "              INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);
1101   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"              [%d] %d \n",lu->myid,lu->id.INFO(23));CHKERRQ(ierr);
1102   ierr = PetscViewerFlush(viewer);
1103 
1104   if (!lu->myid){ /* information from the host */
1105     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr);
1106     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr);
1107     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr);
1108 
1109     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr);
1110     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr);
1111     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr);
1112     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr);
1113     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr);
1114     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr);
1115     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",lu->id.INFOG(9));CHKERRQ(ierr);
1116     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",lu->id.INFOG(10));CHKERRQ(ierr);
1117     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",lu->id.INFOG(11));CHKERRQ(ierr);
1118     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr);
1119     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr);
1120     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr);
1121     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr);
1122     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);
1123     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);
1124     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);
1125     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);
1126      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr);
1127      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);
1128      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);
1129      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",lu->id.INFOG(23));CHKERRQ(ierr);
1130      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",lu->id.INFOG(24));CHKERRQ(ierr);
1131      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",lu->id.INFOG(25));CHKERRQ(ierr);
1132   }
1133   PetscFunctionReturn(0);
1134 }
1135 
1136 #undef __FUNCT__
1137 #define __FUNCT__ "MatView_MUMPS"
1138 PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
1139 {
1140   PetscErrorCode    ierr;
1141   PetscTruth        iascii;
1142   PetscViewerFormat format;
1143   Mat_MUMPS         *mumps=(Mat_MUMPS*)A->spptr;
1144 
1145   PetscFunctionBegin;
1146   /* check if matrix is mumps type */
1147   if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0);
1148 
1149   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1150   if (iascii) {
1151     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1152     if (format == PETSC_VIEWER_ASCII_INFO || mumps->mumpsview){
1153       ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
1154       ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                  %d \n",mumps->id.sym);CHKERRQ(ierr);
1155       ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):           %d \n",mumps->id.par);CHKERRQ(ierr);
1156       if (mumps->mumpsview){ /* View all MUMPS parameters */
1157         ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr);
1158       }
1159     }
1160   }
1161   PetscFunctionReturn(0);
1162 }
1163 
1164 #undef __FUNCT__
1165 #define __FUNCT__ "MatGetInfo_MUMPS"
1166 PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
1167 {
1168   Mat_MUMPS  *mumps =(Mat_MUMPS*)A->spptr;
1169 
1170   PetscFunctionBegin;
1171   info->block_size        = 1.0;
1172   info->nz_allocated      = mumps->id.INFOG(20);
1173   info->nz_used           = mumps->id.INFOG(20);
1174   info->nz_unneeded       = 0.0;
1175   info->assemblies        = 0.0;
1176   info->mallocs           = 0.0;
1177   info->memory            = 0.0;
1178   info->fill_ratio_given  = 0;
1179   info->fill_ratio_needed = 0;
1180   info->factor_mallocs    = 0;
1181   PetscFunctionReturn(0);
1182 }
1183 
1184 /*MC
1185   MATSOLVERMUMPS -  A matrix type providing direct solvers (LU and Cholesky) for
1186   distributed and sequential matrices via the external package MUMPS.
1187 
1188   Works with MATAIJ and MATSBAIJ matrices
1189 
1190   Options Database Keys:
1191 + -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
1192 . -mat_mumps_icntl_4 <0,...,4> - print level
1193 . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
1194 . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
1195 . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
1196 . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
1197 . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
1198 . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
1199 . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
1200 . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
1201 . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
1202 . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
1203 . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
1204 - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
1205 
1206   Level: beginner
1207 
1208 .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage
1209 
1210 M*/
1211 
1212 EXTERN_C_BEGIN
1213 #undef __FUNCT__
1214 #define __FUNCT__ "MatFactorGetSolverPackage_mumps"
1215 PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type)
1216 {
1217   PetscFunctionBegin;
1218   *type = MATSOLVERMUMPS;
1219   PetscFunctionReturn(0);
1220 }
1221 EXTERN_C_END
1222 
1223 EXTERN_C_BEGIN
1224 /* MatGetFactor for Seq and MPI AIJ matrices */
1225 #undef __FUNCT__
1226 #define __FUNCT__ "MatGetFactor_aij_mumps"
1227 PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
1228 {
1229   Mat            B;
1230   PetscErrorCode ierr;
1231   Mat_MUMPS      *mumps;
1232   PetscTruth     isSeqAIJ;
1233 
1234   PetscFunctionBegin;
1235   /* Create the factorization matrix */
1236   ierr = PetscTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr);
1237   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
1238   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
1239   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1240   if (isSeqAIJ) {
1241     ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
1242   } else {
1243     ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
1244   }
1245 
1246   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
1247   B->ops->view             = MatView_MUMPS;
1248   B->ops->getinfo          = MatGetInfo_MUMPS;
1249   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
1250   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl",MatMumpsSetIcntl);CHKERRQ(ierr);
1251   if (ftype == MAT_FACTOR_LU) {
1252     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
1253     B->factortype = MAT_FACTOR_LU;
1254     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
1255     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
1256   } else {
1257     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
1258     B->factortype = MAT_FACTOR_CHOLESKY;
1259     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
1260     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
1261   }
1262 
1263   mumps->CleanUpMUMPS              = PETSC_FALSE;
1264   mumps->isAIJ                     = PETSC_TRUE;
1265   mumps->scat_rhs                  = PETSC_NULL;
1266   mumps->scat_sol                  = PETSC_NULL;
1267   mumps->nSolve                    = 0;
1268   mumps->MatDestroy                = B->ops->destroy;
1269   B->ops->destroy                  = MatDestroy_MUMPS;
1270   B->spptr                         = (void*)mumps;
1271 
1272   *F = B;
1273   PetscFunctionReturn(0);
1274 }
1275 EXTERN_C_END
1276 
1277 
1278 EXTERN_C_BEGIN
1279 /* MatGetFactor for Seq and MPI SBAIJ matrices */
1280 #undef __FUNCT__
1281 #define __FUNCT__ "MatGetFactor_sbaij_mumps"
1282 PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
1283 {
1284   Mat            B;
1285   PetscErrorCode ierr;
1286   Mat_MUMPS      *mumps;
1287   PetscTruth     isSeqSBAIJ;
1288 
1289   PetscFunctionBegin;
1290   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
1291   if(A->rmap->bs > 1) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with block size > 1 with MUMPS Cholesky, use AIJ matrix instead");
1292   ierr = PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
1293   /* Create the factorization matrix */
1294   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
1295   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
1296   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1297   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
1298   if (isSeqSBAIJ) {
1299     ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
1300     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
1301   } else {
1302     ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
1303     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
1304   }
1305 
1306   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
1307   B->ops->view                   = MatView_MUMPS;
1308   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
1309   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl",MatMumpsSetIcntl);CHKERRQ(ierr);
1310   B->factortype                   = MAT_FACTOR_CHOLESKY;
1311 
1312   mumps->CleanUpMUMPS              = PETSC_FALSE;
1313   mumps->isAIJ                     = PETSC_FALSE;
1314   mumps->scat_rhs                  = PETSC_NULL;
1315   mumps->scat_sol                  = PETSC_NULL;
1316   mumps->nSolve                    = 0;
1317   mumps->MatDestroy                = B->ops->destroy;
1318   B->ops->destroy                  = MatDestroy_MUMPS;
1319   B->spptr                         = (void*)mumps;
1320 
1321   *F = B;
1322   PetscFunctionReturn(0);
1323 }
1324 EXTERN_C_END
1325 
1326 EXTERN_C_BEGIN
1327 #undef __FUNCT__
1328 #define __FUNCT__ "MatGetFactor_baij_mumps"
1329 PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
1330 {
1331   Mat            B;
1332   PetscErrorCode ierr;
1333   Mat_MUMPS      *mumps;
1334   PetscTruth     isSeqBAIJ;
1335 
1336   PetscFunctionBegin;
1337   /* Create the factorization matrix */
1338   ierr = PetscTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr);
1339   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
1340   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
1341   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1342   if (isSeqBAIJ) {
1343     ierr = MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,PETSC_NULL);CHKERRQ(ierr);
1344   } else {
1345     ierr = MatMPIBAIJSetPreallocation(B,A->rmap->bs,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
1346   }
1347 
1348   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
1349   if (ftype == MAT_FACTOR_LU) {
1350     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
1351     B->factortype = MAT_FACTOR_LU;
1352     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
1353     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
1354   }
1355   else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n");
1356 
1357   B->ops->view             = MatView_MUMPS;
1358   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
1359   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl",MatMumpsSetIcntl);CHKERRQ(ierr);
1360 
1361   mumps->CleanUpMUMPS              = PETSC_FALSE;
1362   mumps->isAIJ                     = PETSC_TRUE;
1363   mumps->scat_rhs                  = PETSC_NULL;
1364   mumps->scat_sol                  = PETSC_NULL;
1365   mumps->nSolve                    = 0;
1366   mumps->MatDestroy                = B->ops->destroy;
1367   B->ops->destroy                  = MatDestroy_MUMPS;
1368   B->spptr                         = (void*)mumps;
1369 
1370   *F = B;
1371   PetscFunctionReturn(0);
1372 }
1373 EXTERN_C_END
1374 
1375 /* -------------------------------------------------------------------------------------------*/
1376 /*@
1377   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()
1378 
1379    Collective on Mat
1380 
1381    Input Parameters:
1382 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1383 .  idx - index of MUMPS parameter array ICNTL()
1384 -  icntl - value of MUMPS ICNTL(imumps)
1385 
1386   Options Database:
1387 .   -mat_mumps_icntl_<idx> <icntl>
1388 
1389    Level: beginner
1390 
1391    References: MUMPS Users' Guide
1392 
1393 .seealso: MatGetFactor()
1394 @*/
1395 #undef __FUNCT__
1396 #define __FUNCT__ "MatMumpsSetIcntl"
1397 PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt idx,PetscInt icntl)
1398 {
1399   Mat_MUMPS      *lu =(Mat_MUMPS*)(F)->spptr;
1400 
1401   PetscFunctionBegin;
1402   lu->id.ICNTL(idx) = icntl;
1403   PetscFunctionReturn(0);
1404 }
1405 
1406