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