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