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