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