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