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