Actual source code: mumps.c
1: #define PETSCMAT_DLL
3: /*
4: Provides an interface to the MUMPS sparse solver
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
12: #if defined(PETSC_USE_COMPLEX)
13: #include "zmumps_c.h"
14: #else
15: #include "dmumps_c.h"
16: #endif
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 INFO(I) info[(I)-1]
25: #define RINFOG(I) rinfog[(I)-1]
26: #define RINFO(I) rinfo[(I)-1]
28: typedef struct {
29: #if defined(PETSC_USE_COMPLEX)
30: ZMUMPS_STRUC_C id;
31: #else
32: DMUMPS_STRUC_C id;
33: #endif
34: MatStructure matstruc;
35: PetscMPIInt myid,size;
36: PetscInt *irn,*jcn,sym,nSolve;
37: PetscScalar *val;
38: MPI_Comm comm_mumps;
39: VecScatter scat_rhs, scat_sol;
40: PetscTruth isAIJ,CleanUpMUMPS;
41: Vec b_seq,x_seq;
42: PetscErrorCode (*MatDestroy)(Mat);
43: } Mat_MUMPS;
45: EXTERN PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);
47: /* convert Petsc mpiaij matrix to triples: row[nz], col[nz], val[nz] */
48: /*
49: input:
50: A - matrix in mpiaij or mpisbaij (bs=1) format
51: shift - 0: C style output triple; 1: Fortran style output triple.
52: valOnly - FALSE: spaces are allocated and values are set for the triple
53: TRUE: only the values in v array are updated
54: output:
55: nnz - dim of r, c, and v (number of local nonzero entries of A)
56: r, c, v - row and col index, matrix values (matrix triples)
57: */
58: PetscErrorCode MatConvertToTriples(Mat A,int shift,PetscTruth valOnly,int *nnz,int **r, int **c, PetscScalar **v)
59: {
60: PetscInt *ai, *aj, *bi, *bj, rstart,nz, *garray;
62: PetscInt i,j,jj,jB,irow,m=A->rmap->n,*ajj,*bjj,countA,countB,colA_start,jcol;
63: PetscInt *row,*col;
64: PetscScalar *av, *bv,*val;
65: PetscTruth isAIJ;
68: PetscTypeCompare((PetscObject)A,MATMPIAIJ,&isAIJ);
69: if (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= A->rmap->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 (A->rmap->bs > 1) SETERRQ1(PETSC_ERR_SUP," bs=%d is not supported yet\n", A->rmap->bs);
83: nz = aa->nz + bb->nz;
84: ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
85: garray = mat->garray;
86: av=aa->a; bv=bb->a;
87: }
89: if (!valOnly){
90: PetscMalloc(nz*sizeof(PetscInt) ,&row);
91: PetscMalloc(nz*sizeof(PetscInt),&col);
92: PetscMalloc(nz*sizeof(PetscScalar),&val);
93: *r = row; *c = col; *v = val;
94: } else {
95: row = *r; col = *c; val = *v;
96: }
97: *nnz = nz;
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];
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;
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: return(0);
144: }
148: PetscErrorCode MatDestroy_MUMPS(Mat A)
149: {
150: Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr;
152: PetscMPIInt size=lu->size;
155: if (lu->CleanUpMUMPS) {
156: /* Terminate instance, deallocate memories */
157: if (size > 1){
158: PetscFree(lu->id.sol_loc);
159: VecScatterDestroy(lu->scat_rhs);
160: VecDestroy(lu->b_seq);
161: if (lu->nSolve && lu->scat_sol){VecScatterDestroy(lu->scat_sol);}
162: if (lu->nSolve && lu->x_seq){VecDestroy(lu->x_seq);}
163: PetscFree(lu->val);
164: }
165: lu->id.job=JOB_END;
166: #if defined(PETSC_USE_COMPLEX)
167: zmumps_c(&lu->id);
168: #else
169: dmumps_c(&lu->id);
170: #endif
171: PetscFree(lu->irn);
172: PetscFree(lu->jcn);
173: MPI_Comm_free(&(lu->comm_mumps));
174: }
175: (lu->MatDestroy)(A);
176: return(0);
177: }
181: PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x)
182: {
183: Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr;
184: PetscScalar *array;
185: Vec x_seq;
186: IS is_iden,is_petsc;
188: PetscInt i;
191: lu->id.nrhs = 1;
192: x_seq = lu->b_seq;
193: if (lu->size > 1){
194: /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
195: VecScatterBegin(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);
196: VecScatterEnd(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);
197: if (!lu->myid) {VecGetArray(x_seq,&array);}
198: } else { /* size == 1 */
199: VecCopy(b,x);
200: VecGetArray(x,&array);
201: }
202: if (!lu->myid) { /* define rhs on the host */
203: #if defined(PETSC_USE_COMPLEX)
204: lu->id.rhs = (mumps_double_complex*)array;
205: #else
206: lu->id.rhs = array;
207: #endif
208: }
209: if (lu->size == 1){
210: VecRestoreArray(x,&array);
211: } else if (!lu->myid){
212: VecRestoreArray(x_seq,&array);
213: }
215: if (lu->size > 1){
216: /* distributed solution */
217: lu->id.ICNTL(21) = 1;
218: if (!lu->nSolve){
219: /* Create x_seq=sol_loc for repeated use */
220: PetscInt lsol_loc;
221: PetscScalar *sol_loc;
222: lsol_loc = lu->id.INFO(23); /* length of sol_loc */
223: PetscMalloc((1+lsol_loc)*(sizeof(PetscScalar)+sizeof(PetscInt)),&sol_loc);
224: lu->id.isol_loc = (PetscInt *)(sol_loc + lsol_loc);
225: lu->id.lsol_loc = lsol_loc;
226: #if defined(PETSC_USE_COMPLEX)
227: lu->id.sol_loc = (mumps_double_complex*)sol_loc;
228: #else
229: lu->id.sol_loc = sol_loc;
230: #endif
231: VecCreateSeqWithArray(PETSC_COMM_SELF,lsol_loc,sol_loc,&lu->x_seq);
232: }
233: }
235: /* solve phase */
236: /*-------------*/
237: lu->id.job = 3;
238: #if defined(PETSC_USE_COMPLEX)
239: zmumps_c(&lu->id);
240: #else
241: dmumps_c(&lu->id);
242: #endif
243: if (lu->id.INFOG(1) < 0) {
244: SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1));
245: }
247: if (lu->size > 1) { /* convert mumps distributed solution to petsc mpi x */
248: if (!lu->nSolve){ /* create scatter scat_sol */
249: ISCreateStride(PETSC_COMM_SELF,lu->id.lsol_loc,0,1,&is_iden); /* from */
250: for (i=0; i<lu->id.lsol_loc; i++){
251: lu->id.isol_loc[i] -= 1; /* change Fortran style to C style */
252: }
253: ISCreateGeneral(PETSC_COMM_SELF,lu->id.lsol_loc,lu->id.isol_loc,&is_petsc); /* to */
254: VecScatterCreate(lu->x_seq,is_iden,x,is_petsc,&lu->scat_sol);
255: ISDestroy(is_iden);
256: ISDestroy(is_petsc);
257: }
258: VecScatterBegin(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);
259: VecScatterEnd(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);
260: }
261: lu->nSolve++;
262: return(0);
263: }
265: #if !defined(PETSC_USE_COMPLEX)
266: /*
267: input:
268: F: numeric factor
269: output:
270: nneg: total number of negative pivots
271: nzero: 0
272: npos: (global dimension of F) - nneg
273: */
277: PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
278: {
279: Mat_MUMPS *lu =(Mat_MUMPS*)F->spptr;
281: PetscMPIInt size;
284: MPI_Comm_size(((PetscObject)F)->comm,&size);
285: /* MUMPS 4.3.1 calls ScaLAPACK when ICNTL(13)=0 (default), which does not offer the possibility to compute the inertia of a dense matrix. Set ICNTL(13)=1 to skip ScaLAPACK */
286: if (size > 1 && lu->id.ICNTL(13) != 1){
287: SETERRQ1(PETSC_ERR_ARG_WRONG,"ICNTL(13)=%d. -mat_mumps_icntl_13 must be set as 1 for correct global matrix inertia\n",lu->id.INFOG(13));
288: }
289: if (nneg){
290: if (!lu->myid){
291: *nneg = lu->id.INFOG(12);
292: }
293: MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);
294: }
295: if (nzero) *nzero = 0;
296: if (npos) *npos = F->rmap->N - (*nneg);
297: return(0);
298: }
299: #endif /* !defined(PETSC_USE_COMPLEX) */
303: PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
304: {
305: Mat_MUMPS *lu =(Mat_MUMPS*)(F)->spptr;
307: PetscInt rnz,nnz,nz=0,i,M=A->rmap->N,*ai,*aj,icntl;
308: PetscTruth valOnly,flg;
309: Mat F_diag;
310: IS is_iden;
311: Vec b;
312: PetscTruth isSeqAIJ,isSeqSBAIJ;
315: PetscTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);
316: PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);
317: if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
318: (F)->ops->solve = MatSolve_MUMPS;
320: /* Initialize a MUMPS instance */
321: MPI_Comm_rank(((PetscObject)A)->comm, &lu->myid);
322: MPI_Comm_size(((PetscObject)A)->comm,&lu->size);
323: lu->id.job = JOB_INIT;
324: MPI_Comm_dup(((PetscObject)A)->comm,&(lu->comm_mumps));
325: lu->id.comm_fortran = MPI_Comm_c2f(lu->comm_mumps);
327: /* Set mumps options */
328: PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"MUMPS Options","Mat");
329: lu->id.par=1; /* host participates factorizaton and solve */
330: lu->id.sym=lu->sym;
331: if (lu->sym == 2){
332: PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);
333: if (flg && icntl == 1) lu->id.sym=icntl; /* matrix is spd */
334: }
335: #if defined(PETSC_USE_COMPLEX)
336: zmumps_c(&lu->id);
337: #else
338: dmumps_c(&lu->id);
339: #endif
340:
341: if (isSeqAIJ || isSeqSBAIJ){
342: lu->id.ICNTL(18) = 0; /* centralized assembled matrix input */
343: } else {
344: lu->id.ICNTL(18) = 3; /* distributed assembled matrix input */
345: }
347: icntl=-1;
348: lu->id.ICNTL(4) = 0; /* level of printing; overwrite mumps default ICNTL(4)=2 */
349: PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);
350: if ((flg && icntl > 0) || PetscLogPrintInfo) {
351: lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */
352: } else { /* no output */
353: lu->id.ICNTL(1) = 0; /* error message, default= 6 */
354: lu->id.ICNTL(2) = 0; /* output stream for diagnostic printing, statistics, and warning. default=0 */
355: lu->id.ICNTL(3) = 0; /* output stream for global information, default=6 */
356: }
357: PetscOptionsInt("-mat_mumps_icntl_6","ICNTL(6): column permutation and/or scaling to get a zero-free diagonal (0 to 7)","None",lu->id.ICNTL(6),&lu->id.ICNTL(6),PETSC_NULL);
358: icntl=-1;
359: PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);
360: if (flg) {
361: if (icntl== 1){
362: SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n");
363: } else {
364: lu->id.ICNTL(7) = icntl;
365: }
366: }
367: PetscOptionsInt("-mat_mumps_icntl_8","ICNTL(8): scaling strategy (-2 to 7 or 77)","None",lu->id.ICNTL(8),&lu->id.ICNTL(8),PETSC_NULL);
368: 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);
369: PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",lu->id.ICNTL(10),&lu->id.ICNTL(10),PETSC_NULL);
370: PetscOptionsInt("-mat_mumps_icntl_11","ICNTL(11): statistics related to the linear system solved (via -ksp_view)","None",lu->id.ICNTL(11),&lu->id.ICNTL(11),PETSC_NULL);
371: PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control: defines the ordering strategy with scaling constraints (0 to 3","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);
372: PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control: with or without ScaLAPACK","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);
373: PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage of estimated workspace increase","None",lu->id.ICNTL(14),&lu->id.ICNTL(14),PETSC_NULL);
374: PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): Schur complement","None",lu->id.ICNTL(19),&lu->id.ICNTL(19),PETSC_NULL);
376: PetscOptionsInt("-mat_mumps_icntl_22","ICNTL(22): in-core/out-of-core facility (0 or 1)","None",lu->id.ICNTL(22),&lu->id.ICNTL(22),PETSC_NULL);
377: PetscOptionsInt("-mat_mumps_icntl_23","ICNTL(23): max size of the working memory (MB) that can allocate per processor","None",lu->id.ICNTL(23),&lu->id.ICNTL(23),PETSC_NULL);
378: PetscOptionsInt("-mat_mumps_icntl_24","ICNTL(24): detection of null pivot rows (0 or 1)","None",lu->id.ICNTL(24),&lu->id.ICNTL(24),PETSC_NULL);
379: PetscOptionsInt("-mat_mumps_icntl_25","ICNTL(25): computation of a null space basis","None",lu->id.ICNTL(25),&lu->id.ICNTL(25),PETSC_NULL);
380: PetscOptionsInt("-mat_mumps_icntl_26","ICNTL(26): Schur options for right-hand side or solution vector","None",lu->id.ICNTL(26),&lu->id.ICNTL(26),PETSC_NULL);
381: PetscOptionsInt("-mat_mumps_icntl_27","ICNTL(27): experimental parameter","None",lu->id.ICNTL(27),&lu->id.ICNTL(27),PETSC_NULL);
383: PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",lu->id.CNTL(1),&lu->id.CNTL(1),PETSC_NULL);
384: PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",lu->id.CNTL(2),&lu->id.CNTL(2),PETSC_NULL);
385: PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",lu->id.CNTL(3),&lu->id.CNTL(3),PETSC_NULL);
386: PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",lu->id.CNTL(4),&lu->id.CNTL(4),PETSC_NULL);
387: PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",lu->id.CNTL(5),&lu->id.CNTL(5),PETSC_NULL);
388: PetscOptionsEnd();
389: }
391: /* define matrix A */
392: switch (lu->id.ICNTL(18)){
393: case 0: /* centralized assembled matrix input (size=1) */
394: if (!lu->myid) {
395: if (isSeqAIJ){
396: Mat_SeqAIJ *aa = (Mat_SeqAIJ*)A->data;
397: nz = aa->nz;
398: ai = aa->i; aj = aa->j; lu->val = aa->a;
399: } else if (isSeqSBAIJ) {
400: Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data;
401: nz = aa->nz;
402: ai = aa->i; aj = aa->j; lu->val = aa->a;
403: } else {
404: SETERRQ(PETSC_ERR_SUP,"No mumps factorization for this matrix type");
405: }
406: if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */
407: PetscMalloc(nz*sizeof(PetscInt),&lu->irn);
408: PetscMalloc(nz*sizeof(PetscInt),&lu->jcn);
409: nz = 0;
410: for (i=0; i<M; i++){
411: rnz = ai[i+1] - ai[i];
412: while (rnz--) { /* Fortran row/col index! */
413: lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++;
414: }
415: }
416: }
417: }
418: break;
419: case 3: /* distributed assembled matrix input (size>1) */
420: if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
421: valOnly = PETSC_FALSE;
422: } else {
423: valOnly = PETSC_TRUE; /* only update mat values, not row and col index */
424: }
425: MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);
426: break;
427: default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS.");
428: }
430: /* analysis phase */
431: /*----------------*/
432: if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
433: lu->id.job = 1;
435: lu->id.n = M;
436: switch (lu->id.ICNTL(18)){
437: case 0: /* centralized assembled matrix input */
438: if (!lu->myid) {
439: lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
440: if (lu->id.ICNTL(6)>1){
441: #if defined(PETSC_USE_COMPLEX)
442: lu->id.a = (mumps_double_complex*)lu->val;
443: #else
444: lu->id.a = lu->val;
445: #endif
446: }
447: }
448: break;
449: case 3: /* distributed assembled matrix input (size>1) */
450: lu->id.nz_loc = nnz;
451: lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
452: if (lu->id.ICNTL(6)>1) {
453: #if defined(PETSC_USE_COMPLEX)
454: lu->id.a_loc = (mumps_double_complex*)lu->val;
455: #else
456: lu->id.a_loc = lu->val;
457: #endif
458: }
459: /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
460: if (!lu->myid){
461: VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);
462: ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);
463: } else {
464: VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);
465: ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);
466: }
467: VecCreate(((PetscObject)A)->comm,&b);
468: VecSetSizes(b,A->rmap->n,PETSC_DECIDE);
469: VecSetFromOptions(b);
471: VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);
472: ISDestroy(is_iden);
473: VecDestroy(b);
474: break;
475: }
476: #if defined(PETSC_USE_COMPLEX)
477: zmumps_c(&lu->id);
478: #else
479: dmumps_c(&lu->id);
480: #endif
481: if (lu->id.INFOG(1) < 0) {
482: SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
483: }
484: }
486: /* numerical factorization phase */
487: /*-------------------------------*/
488: lu->id.job = 2;
489: if(!lu->id.ICNTL(18)) {
490: if (!lu->myid) {
491: #if defined(PETSC_USE_COMPLEX)
492: lu->id.a = (mumps_double_complex*)lu->val;
493: #else
494: lu->id.a = lu->val;
495: #endif
496: }
497: } else {
498: #if defined(PETSC_USE_COMPLEX)
499: lu->id.a_loc = (mumps_double_complex*)lu->val;
500: #else
501: lu->id.a_loc = lu->val;
502: #endif
503: }
504: #if defined(PETSC_USE_COMPLEX)
505: zmumps_c(&lu->id);
506: #else
507: dmumps_c(&lu->id);
508: #endif
509: if (lu->id.INFOG(1) < 0) {
510: if (lu->id.INFO(1) == -13) {
511: SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",lu->id.INFO(2));
512: } else {
513: SETERRQ2(PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: INFO(1)=%d, INFO(2)=%d\n",lu->id.INFO(1),lu->id.INFO(2));
514: }
515: }
517: if (!lu->myid && lu->id.ICNTL(16) > 0){
518: SETERRQ1(PETSC_ERR_LIB," lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16));
519: }
521: if (lu->size > 1){
522: if ((F)->factor == MAT_FACTOR_LU){
523: F_diag = ((Mat_MPIAIJ *)(F)->data)->A;
524: } else {
525: F_diag = ((Mat_MPISBAIJ *)(F)->data)->A;
526: }
527: F_diag->assembled = PETSC_TRUE;
528: if (lu->nSolve){
529: VecScatterDestroy(lu->scat_sol);
530: PetscFree(lu->id.sol_loc);
531: VecDestroy(lu->x_seq);
532: }
533: }
534: (F)->assembled = PETSC_TRUE;
535: lu->matstruc = SAME_NONZERO_PATTERN;
536: lu->CleanUpMUMPS = PETSC_TRUE;
537: lu->nSolve = 0;
538: return(0);
539: }
542: /* Note the Petsc r and c permutations are ignored */
545: PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
546: {
547: Mat_MUMPS *lu = (Mat_MUMPS*)F->spptr;
550: lu->sym = 0;
551: lu->matstruc = DIFFERENT_NONZERO_PATTERN;
552: F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
553: return(0);
554: }
557: /* Note the Petsc r permutation is ignored */
560: PetscErrorCode MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
561: {
562: Mat_MUMPS *lu = (Mat_MUMPS*)(F)->spptr;
565: lu->sym = 2;
566: lu->matstruc = DIFFERENT_NONZERO_PATTERN;
567: (F)->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
568: #if !defined(PETSC_USE_COMPLEX)
569: (F)->ops->getinertia = MatGetInertia_SBAIJMUMPS;
570: #endif
571: return(0);
572: }
576: PetscErrorCode MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) {
577: Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr;
581: /* check if matrix is mumps type */
582: if (A->ops->solve != MatSolve_MUMPS) return(0);
584: PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");
585: PetscViewerASCIIPrintf(viewer," SYM (matrix type): %d \n",lu->id.sym);
586: PetscViewerASCIIPrintf(viewer," PAR (host participation): %d \n",lu->id.par);
587: PetscViewerASCIIPrintf(viewer," ICNTL(1) (output for error): %d \n",lu->id.ICNTL(1));
588: PetscViewerASCIIPrintf(viewer," ICNTL(2) (output of diagnostic msg):%d \n",lu->id.ICNTL(2));
589: PetscViewerASCIIPrintf(viewer," ICNTL(3) (output for global info): %d \n",lu->id.ICNTL(3));
590: PetscViewerASCIIPrintf(viewer," ICNTL(4) (level of printing): %d \n",lu->id.ICNTL(4));
591: PetscViewerASCIIPrintf(viewer," ICNTL(5) (input mat struct): %d \n",lu->id.ICNTL(5));
592: PetscViewerASCIIPrintf(viewer," ICNTL(6) (matrix prescaling): %d \n",lu->id.ICNTL(6));
593: PetscViewerASCIIPrintf(viewer," ICNTL(7) (matrix ordering): %d \n",lu->id.ICNTL(7));
594: PetscViewerASCIIPrintf(viewer," ICNTL(8) (scalling strategy): %d \n",lu->id.ICNTL(8));
595: PetscViewerASCIIPrintf(viewer," ICNTL(9) (A/A^T x=b is solved): %d \n",lu->id.ICNTL(9));
596: PetscViewerASCIIPrintf(viewer," ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));
597: PetscViewerASCIIPrintf(viewer," ICNTL(11) (error analysis): %d \n",lu->id.ICNTL(11));
598: if (!lu->myid && lu->id.ICNTL(11)>0) {
599: PetscPrintf(PETSC_COMM_SELF," RINFOG(4) (inf norm of input mat): %g\n",lu->id.RINFOG(4));
600: PetscPrintf(PETSC_COMM_SELF," RINFOG(5) (inf norm of solution): %g\n",lu->id.RINFOG(5));
601: PetscPrintf(PETSC_COMM_SELF," RINFOG(6) (inf norm of residual): %g\n",lu->id.RINFOG(6));
602: PetscPrintf(PETSC_COMM_SELF," RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",lu->id.RINFOG(7),lu->id.RINFOG(8));
603: PetscPrintf(PETSC_COMM_SELF," RINFOG(9) (error estimate): %g \n",lu->id.RINFOG(9));
604: PetscPrintf(PETSC_COMM_SELF," RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",lu->id.RINFOG(10),lu->id.RINFOG(11));
605:
606: }
607: PetscViewerASCIIPrintf(viewer," ICNTL(12) (efficiency control): %d \n",lu->id.ICNTL(12));
608: PetscViewerASCIIPrintf(viewer," ICNTL(13) (efficiency control): %d \n",lu->id.ICNTL(13));
609: PetscViewerASCIIPrintf(viewer," ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));
610: /* ICNTL(15-17) not used */
611: PetscViewerASCIIPrintf(viewer," ICNTL(18) (input mat struct): %d \n",lu->id.ICNTL(18));
612: PetscViewerASCIIPrintf(viewer," ICNTL(19) (Shur complement info): %d \n",lu->id.ICNTL(19));
613: PetscViewerASCIIPrintf(viewer," ICNTL(20) (rhs sparse pattern): %d \n",lu->id.ICNTL(20));
614: PetscViewerASCIIPrintf(viewer," ICNTL(21) (solution struct): %d \n",lu->id.ICNTL(21));
615: PetscViewerASCIIPrintf(viewer," ICNTL(22) (in-core/out-of-core facility): %d \n",lu->id.ICNTL(22));
616: PetscViewerASCIIPrintf(viewer," ICNTL(23) (max size of memory can be allocated locally):%d \n",lu->id.ICNTL(23));
618: PetscViewerASCIIPrintf(viewer," ICNTL(24) (detection of null pivot rows): %d \n",lu->id.ICNTL(24));
619: PetscViewerASCIIPrintf(viewer," ICNTL(25) (computation of a null space basis): %d \n",lu->id.ICNTL(25));
620: PetscViewerASCIIPrintf(viewer," ICNTL(26) (Schur options for rhs or solution): %d \n",lu->id.ICNTL(26));
621: PetscViewerASCIIPrintf(viewer," ICNTL(27) (experimental parameter): %d \n",lu->id.ICNTL(27));
623: PetscViewerASCIIPrintf(viewer," CNTL(1) (relative pivoting threshold): %g \n",lu->id.CNTL(1));
624: PetscViewerASCIIPrintf(viewer," CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));
625: PetscViewerASCIIPrintf(viewer," CNTL(3) (absolute pivoting threshold): %g \n",lu->id.CNTL(3));
626: PetscViewerASCIIPrintf(viewer," CNTL(4) (value of static pivoting): %g \n",lu->id.CNTL(4));
627: PetscViewerASCIIPrintf(viewer," CNTL(5) (fixation for null pivots): %g \n",lu->id.CNTL(5));
629: /* infomation local to each processor */
630: if (!lu->myid) {PetscPrintf(PETSC_COMM_SELF, " RINFO(1) (local estimated flops for the elimination after analysis): \n");}
631: PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %g \n",lu->myid,lu->id.RINFO(1));
632: PetscSynchronizedFlush(((PetscObject)A)->comm);
633: if (!lu->myid) {PetscPrintf(PETSC_COMM_SELF, " RINFO(2) (local estimated flops for the assembly after factorization): \n");}
634: PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %g \n",lu->myid,lu->id.RINFO(2));
635: PetscSynchronizedFlush(((PetscObject)A)->comm);
636: if (!lu->myid) {PetscPrintf(PETSC_COMM_SELF, " RINFO(3) (local estimated flops for the elimination after factorization): \n");}
637: PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %g \n",lu->myid,lu->id.RINFO(3));
638: PetscSynchronizedFlush(((PetscObject)A)->comm);
640: if (!lu->myid) {PetscPrintf(PETSC_COMM_SELF, " INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");}
641: PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %d \n",lu->myid,lu->id.INFO(15));
642: PetscSynchronizedFlush(((PetscObject)A)->comm);
644: if (!lu->myid) {PetscPrintf(PETSC_COMM_SELF, " INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");}
645: PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %d \n",lu->myid,lu->id.INFO(16));
646: PetscSynchronizedFlush(((PetscObject)A)->comm);
648: if (!lu->myid) {PetscPrintf(PETSC_COMM_SELF, " INFO(23) (num of pivots eliminated on this processor after factorization): \n");}
649: PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %d \n",lu->myid,lu->id.INFO(23));
650: PetscSynchronizedFlush(((PetscObject)A)->comm);
652: if (!lu->myid){ /* information from the host */
653: PetscViewerASCIIPrintf(viewer," RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));
654: PetscViewerASCIIPrintf(viewer," RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));
655: PetscViewerASCIIPrintf(viewer," RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));
657: PetscViewerASCIIPrintf(viewer," INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));
658: PetscViewerASCIIPrintf(viewer," INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));
659: PetscViewerASCIIPrintf(viewer," INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));
660: PetscViewerASCIIPrintf(viewer," INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));
661: PetscViewerASCIIPrintf(viewer," INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));
662: PetscViewerASCIIPrintf(viewer," INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));
663: PetscViewerASCIIPrintf(viewer," INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",lu->id.INFOG(9));
664: PetscViewerASCIIPrintf(viewer," INFOG(10) (total integer space store the matrix factors after factorization): %d \n",lu->id.INFOG(10));
665: PetscViewerASCIIPrintf(viewer," INFOG(11) (order of largest frontal matrix after factorization): %d \n",lu->id.INFOG(11));
666: PetscViewerASCIIPrintf(viewer," INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));
667: PetscViewerASCIIPrintf(viewer," INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));
668: PetscViewerASCIIPrintf(viewer," INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));
669: PetscViewerASCIIPrintf(viewer," INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));
670: PetscViewerASCIIPrintf(viewer," INFOG(16) (estimated size (in MB) of all MUMPS internal data for factorization after analysis: value on the most memory consuming processor): %d \n",lu->id.INFOG(16));
671: 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));
672: 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));
673: PetscViewerASCIIPrintf(viewer," INFOG(19) (size of all MUMPS internal data allocated during factorization: sum over all processors): %d \n",lu->id.INFOG(19));
674: PetscViewerASCIIPrintf(viewer," INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));
675: PetscViewerASCIIPrintf(viewer," INFOG(21) (size in MB of memory effectively used during factorization - value on the most memory consuming processor): %d \n",lu->id.INFOG(21));
676: PetscViewerASCIIPrintf(viewer," INFOG(22) (size in MB of memory effectively used during factorization - sum over all processors): %d \n",lu->id.INFOG(22));
677: PetscViewerASCIIPrintf(viewer," INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",lu->id.INFOG(23));
678: PetscViewerASCIIPrintf(viewer," INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",lu->id.INFOG(24));
679: PetscViewerASCIIPrintf(viewer," INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",lu->id.INFOG(25));
680: }
682: return(0);
683: }
687: PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
688: {
689: PetscErrorCode ierr;
690: PetscTruth iascii;
691: PetscViewerFormat format;
694: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
695: if (iascii) {
696: PetscViewerGetFormat(viewer,&format);
697: if (format == PETSC_VIEWER_ASCII_INFO){
698: MatFactorInfo_MUMPS(A,viewer);
699: }
700: }
701: return(0);
702: }
706: PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
707: {
708: Mat_MUMPS *lu =(Mat_MUMPS*)A->spptr;
711: info->block_size = 1.0;
712: info->nz_allocated = lu->id.INFOG(20);
713: info->nz_used = lu->id.INFOG(20);
714: info->nz_unneeded = 0.0;
715: info->assemblies = 0.0;
716: info->mallocs = 0.0;
717: info->memory = 0.0;
718: info->fill_ratio_given = 0;
719: info->fill_ratio_needed = 0;
720: info->factor_mallocs = 0;
721: return(0);
722: }
724: /*MC
725: MAT_SOLVER_MUMPS - A matrix type providing direct solvers (LU and Cholesky) for
726: distributed and sequential matrices via the external package MUMPS.
728: Works with MATAIJ and MATSBAIJ matrices
730: Options Database Keys:
731: + -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
732: . -mat_mumps_icntl_4 <0,...,4> - print level
733: . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
734: . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
735: . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
736: . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
737: . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
738: . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
739: . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
740: . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
741: . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
742: . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
743: . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
744: - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
746: Level: beginner
748: .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage
750: M*/
755: PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type)
756: {
758: *type = MAT_SOLVER_MUMPS;
759: return(0);
760: }
764: /*
765: The seq and mpi versions of this function are the same
766: */
769: PetscErrorCode MatGetFactor_seqaij_mumps(Mat A,MatFactorType ftype,Mat *F)
770: {
771: Mat B;
773: Mat_MUMPS *mumps;
776: if (ftype != MAT_FACTOR_LU) {
777: SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with MUMPS Cholesky, use SBAIJ matrix");
778: }
779: /* Create the factorization matrix */
780: MatCreate(((PetscObject)A)->comm,&B);
781: MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
782: MatSetType(B,((PetscObject)A)->type_name);
783: MatSeqAIJSetPreallocation(B,0,PETSC_NULL);
785: B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
786: B->ops->view = MatView_MUMPS;
787: B->ops->getinfo = MatGetInfo_MUMPS;
788: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);
789: B->factor = MAT_FACTOR_LU;
791: PetscNewLog(B,Mat_MUMPS,&mumps);
792: mumps->CleanUpMUMPS = PETSC_FALSE;
793: mumps->isAIJ = PETSC_TRUE;
794: mumps->scat_rhs = PETSC_NULL;
795: mumps->scat_sol = PETSC_NULL;
796: mumps->nSolve = 0;
797: mumps->MatDestroy = B->ops->destroy;
798: B->ops->destroy = MatDestroy_MUMPS;
799: B->spptr = (void*)mumps;
801: *F = B;
802: return(0);
803: }
809: PetscErrorCode MatGetFactor_mpiaij_mumps(Mat A,MatFactorType ftype,Mat *F)
810: {
811: Mat B;
813: Mat_MUMPS *mumps;
816: if (ftype != MAT_FACTOR_LU) {
817: SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with MUMPS Cholesky, use SBAIJ matrix");
818: }
819: /* Create the factorization matrix */
820: MatCreate(((PetscObject)A)->comm,&B);
821: MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
822: MatSetType(B,((PetscObject)A)->type_name);
823: MatSeqAIJSetPreallocation(B,0,PETSC_NULL);
824: MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);
826: B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
827: B->ops->view = MatView_MUMPS;
828: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);
829: B->factor = MAT_FACTOR_LU;
831: PetscNewLog(B,Mat_MUMPS,&mumps);
832: mumps->CleanUpMUMPS = PETSC_FALSE;
833: mumps->isAIJ = PETSC_TRUE;
834: mumps->scat_rhs = PETSC_NULL;
835: mumps->scat_sol = PETSC_NULL;
836: mumps->nSolve = 0;
837: mumps->MatDestroy = B->ops->destroy;
838: B->ops->destroy = MatDestroy_MUMPS;
839: B->spptr = (void*)mumps;
841: *F = B;
842: return(0);
843: }
849: PetscErrorCode MatGetFactor_seqsbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
850: {
851: Mat B;
853: Mat_MUMPS *mumps;
856: if (ftype != MAT_FACTOR_CHOLESKY) {
857: SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
858: }
859: /* Create the factorization matrix */
860: MatCreate(((PetscObject)A)->comm,&B);
861: MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
862: MatSetType(B,((PetscObject)A)->type_name);
863: MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);
864: MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);
866: B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS;
867: B->ops->view = MatView_MUMPS;
868: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);
870: B->factor = MAT_FACTOR_CHOLESKY;
872: PetscNewLog(B,Mat_MUMPS,&mumps);
873: mumps->CleanUpMUMPS = PETSC_FALSE;
874: mumps->isAIJ = PETSC_TRUE;
875: mumps->scat_rhs = PETSC_NULL;
876: mumps->scat_sol = PETSC_NULL;
877: mumps->nSolve = 0;
878: mumps->MatDestroy = B->ops->destroy;
879: B->ops->destroy = MatDestroy_MUMPS;
880: B->spptr = (void*)mumps;
882: *F = B;
883: return(0);
884: }
890: PetscErrorCode MatGetFactor_mpisbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
891: {
892: Mat B;
894: Mat_MUMPS *mumps;
897: if (ftype != MAT_FACTOR_CHOLESKY) {
898: SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
899: }
900: /* Create the factorization matrix */
901: MatCreate(((PetscObject)A)->comm,&B);
902: MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
903: MatSetType(B,((PetscObject)A)->type_name);
904: MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);
905: MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);
907: B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS;
908: B->ops->view = MatView_MUMPS;
909: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);
910: B->factor = MAT_FACTOR_CHOLESKY;
912: PetscNewLog(B,Mat_MUMPS,&mumps);
913: mumps->CleanUpMUMPS = PETSC_FALSE;
914: mumps->isAIJ = PETSC_TRUE;
915: mumps->scat_rhs = PETSC_NULL;
916: mumps->scat_sol = PETSC_NULL;
917: mumps->nSolve = 0;
918: mumps->MatDestroy = B->ops->destroy;
919: B->ops->destroy = MatDestroy_MUMPS;
920: B->spptr = (void*)mumps;
922: *F = B;
923: return(0);
924: }