Actual source code: mumps.c

petsc-master 2018-02-22
Report Typos and Errors

  2: /*
  3:     Provides an interface to the MUMPS sparse solver
  4: */

  6:  #include <../src/mat/impls/aij/mpi/mpiaij.h>
  7:  #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>

  9: EXTERN_C_BEGIN
 10: #if defined(PETSC_USE_COMPLEX)
 11: #if defined(PETSC_USE_REAL_SINGLE)
 12: #include <cmumps_c.h>
 13: #else
 14: #include <zmumps_c.h>
 15: #endif
 16: #else
 17: #if defined(PETSC_USE_REAL_SINGLE)
 18: #include <smumps_c.h>
 19: #else
 20: #include <dmumps_c.h>
 21: #endif
 22: #endif
 23: EXTERN_C_END
 24: #define JOB_INIT -1
 25: #define JOB_FACTSYMBOLIC 1
 26: #define JOB_FACTNUMERIC 2
 27: #define JOB_SOLVE 3
 28: #define JOB_END -2

 30: /* calls to MUMPS */
 31: #if defined(PETSC_USE_COMPLEX)
 32: #if defined(PETSC_USE_REAL_SINGLE)
 33: #define PetscMUMPS_c cmumps_c
 34: #else
 35: #define PetscMUMPS_c zmumps_c
 36: #endif
 37: #else
 38: #if defined(PETSC_USE_REAL_SINGLE)
 39: #define PetscMUMPS_c smumps_c
 40: #else
 41: #define PetscMUMPS_c dmumps_c
 42: #endif
 43: #endif

 45: /* declare MumpsScalar */
 46: #if defined(PETSC_USE_COMPLEX)
 47: #if defined(PETSC_USE_REAL_SINGLE)
 48: #define MumpsScalar mumps_complex
 49: #else
 50: #define MumpsScalar mumps_double_complex
 51: #endif
 52: #else
 53: #define MumpsScalar PetscScalar
 54: #endif

 56: /* macros s.t. indices match MUMPS documentation */
 57: #define ICNTL(I) icntl[(I)-1]
 58: #define CNTL(I) cntl[(I)-1]
 59: #define INFOG(I) infog[(I)-1]
 60: #define INFO(I) info[(I)-1]
 61: #define RINFOG(I) rinfog[(I)-1]
 62: #define RINFO(I) rinfo[(I)-1]

 64: typedef struct {
 65: #if defined(PETSC_USE_COMPLEX)
 66: #if defined(PETSC_USE_REAL_SINGLE)
 67:   CMUMPS_STRUC_C id;
 68: #else
 69:   ZMUMPS_STRUC_C id;
 70: #endif
 71: #else
 72: #if defined(PETSC_USE_REAL_SINGLE)
 73:   SMUMPS_STRUC_C id;
 74: #else
 75:   DMUMPS_STRUC_C id;
 76: #endif
 77: #endif

 79:   MatStructure matstruc;
 80:   PetscMPIInt  myid,size;
 81:   PetscInt     *irn,*jcn,nz,sym;
 82:   PetscScalar  *val;
 83:   MPI_Comm     comm_mumps;
 84:   PetscBool    isAIJ;
 85:   PetscInt     ICNTL9_pre;           /* check if ICNTL(9) is changed from previous MatSolve */
 86:   VecScatter   scat_rhs, scat_sol;   /* used by MatSolve() */
 87:   Vec          b_seq,x_seq;
 88:   PetscInt     ninfo,*info;          /* display INFO */
 89:   PetscInt     sizeredrhs;
 90:   PetscScalar  *schur_sol;
 91:   PetscInt     schur_sizesol;

 93:   PetscErrorCode (*ConvertToTriples)(Mat, int, MatReuse, int*, int**, int**, PetscScalar**);
 94: } Mat_MUMPS;

 96: extern PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);

 98: static PetscErrorCode MatMumpsResetSchur_Private(Mat_MUMPS* mumps)
 99: {

103:   PetscFree2(mumps->id.listvar_schur,mumps->id.schur);
104:   PetscFree(mumps->id.redrhs);
105:   PetscFree(mumps->schur_sol);
106:   mumps->id.size_schur = 0;
107:   mumps->id.schur_lld  = 0;
108:   mumps->id.ICNTL(19)  = 0;
109:   return(0);
110: }

112: /* solve with rhs in mumps->id.redrhs and return in the same location */
113: static PetscErrorCode MatMumpsSolveSchur_Private(Mat F)
114: {
115:   Mat_MUMPS            *mumps=(Mat_MUMPS*)F->data;
116:   Mat                  S,B,X;
117:   MatFactorSchurStatus schurstatus;
118:   PetscInt             sizesol;
119:   PetscErrorCode       ierr;

122:   MatFactorFactorizeSchurComplement(F);
123:   MatFactorGetSchurComplement(F,&S,&schurstatus);
124:   MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.nrhs,(PetscScalar*)mumps->id.redrhs,&B);
125:   switch (schurstatus) {
126:   case MAT_FACTOR_SCHUR_FACTORED:
127:     MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.nrhs,(PetscScalar*)mumps->id.redrhs,&X);
128:     if (!mumps->id.ICNTL(9)) { /* transpose solve */
129:       MatMatSolveTranspose(S,B,X);
130:     } else {
131:       MatMatSolve(S,B,X);
132:     }
133:     break;
134:   case MAT_FACTOR_SCHUR_INVERTED:
135:     sizesol = mumps->id.nrhs*mumps->id.size_schur;
136:     if (!mumps->schur_sol || sizesol > mumps->schur_sizesol) {
137:       PetscFree(mumps->schur_sol);
138:       PetscMalloc1(sizesol,&mumps->schur_sol);
139:       mumps->schur_sizesol = sizesol;
140:     }
141:     MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.nrhs,mumps->schur_sol,&X);
142:     if (!mumps->id.ICNTL(9)) { /* transpose solve */
143:       MatTransposeMatMult(S,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);
144:     } else {
145:       MatMatMult(S,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);
146:     }
147:     MatCopy(X,B,SAME_NONZERO_PATTERN);
148:     break;
149:   default:
150:     SETERRQ1(PetscObjectComm((PetscObject)F),PETSC_ERR_SUP,"Unhandled MatFactorSchurStatus %D",F->schur_status);
151:     break;
152:   }
153:   MatFactorRestoreSchurComplement(F,&S,schurstatus);
154:   MatDestroy(&B);
155:   MatDestroy(&X);
156:   return(0);
157: }

159: static PetscErrorCode MatMumpsHandleSchur_Private(Mat F, PetscBool expansion)
160: {
161:   Mat_MUMPS     *mumps=(Mat_MUMPS*)F->data;

165:   if (!mumps->id.ICNTL(19)) { /* do nothing when Schur complement has not been computed */
166:     return(0);
167:   }
168:   if (!expansion) { /* prepare for the condensation step */
169:     PetscInt sizeredrhs = mumps->id.nrhs*mumps->id.size_schur;
170:     /* allocate MUMPS internal array to store reduced right-hand sides */
171:     if (!mumps->id.redrhs || sizeredrhs > mumps->sizeredrhs) {
172:       PetscFree(mumps->id.redrhs);
173:       mumps->id.lredrhs = mumps->id.size_schur;
174:       PetscMalloc1(mumps->id.nrhs*mumps->id.lredrhs,&mumps->id.redrhs);
175:       mumps->sizeredrhs = mumps->id.nrhs*mumps->id.lredrhs;
176:     }
177:     mumps->id.ICNTL(26) = 1; /* condensation phase */
178:   } else { /* prepare for the expansion step */
179:     /* solve Schur complement (this has to be done by the MUMPS user, so basically us) */
180:     MatMumpsSolveSchur_Private(F);
181:     mumps->id.ICNTL(26) = 2; /* expansion phase */
182:     PetscMUMPS_c(&mumps->id);
183:     if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
184:     /* restore defaults */
185:     mumps->id.ICNTL(26) = -1;
186:     /* free MUMPS internal array for redrhs if we have solved for multiple rhs in order to save memory space */
187:     if (mumps->id.nrhs > 1) {
188:       PetscFree(mumps->id.redrhs);
189:       mumps->id.lredrhs = 0;
190:       mumps->sizeredrhs = 0;
191:     }
192:   }
193:   return(0);
194: }

196: /*
197:   MatConvertToTriples_A_B - convert Petsc matrix to triples: row[nz], col[nz], val[nz] 

199:   input:
200:     A       - matrix in aij,baij or sbaij (bs=1) format
201:     shift   - 0: C style output triple; 1: Fortran style output triple.
202:     reuse   - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple
203:               MAT_REUSE_MATRIX:   only the values in v array are updated
204:   output:
205:     nnz     - dim of r, c, and v (number of local nonzero entries of A)
206:     r, c, v - row and col index, matrix values (matrix triples)

208:   The returned values r, c, and sometimes v are obtained in a single PetscMalloc(). Then in MatDestroy_MUMPS() it is
209:   freed with PetscFree((mumps->irn);  This is not ideal code, the fact that v is ONLY sometimes part of mumps->irn means
210:   that the PetscMalloc() cannot easily be replaced with a PetscMalloc3(). 

212:  */

214: PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
215: {
216:   const PetscInt *ai,*aj,*ajj,M=A->rmap->n;
217:   PetscInt       nz,rnz,i,j;
219:   PetscInt       *row,*col;
220:   Mat_SeqAIJ     *aa=(Mat_SeqAIJ*)A->data;

223:   *v=aa->a;
224:   if (reuse == MAT_INITIAL_MATRIX) {
225:     nz   = aa->nz;
226:     ai   = aa->i;
227:     aj   = aa->j;
228:     *nnz = nz;
229:     PetscMalloc1(2*nz, &row);
230:     col  = row + nz;

232:     nz = 0;
233:     for (i=0; i<M; i++) {
234:       rnz = ai[i+1] - ai[i];
235:       ajj = aj + ai[i];
236:       for (j=0; j<rnz; j++) {
237:         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
238:       }
239:     }
240:     *r = row; *c = col;
241:   }
242:   return(0);
243: }

245: PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
246: {
247:   Mat_SeqBAIJ    *aa=(Mat_SeqBAIJ*)A->data;
248:   const PetscInt *ai,*aj,*ajj,bs2 = aa->bs2;
249:   PetscInt       bs,M,nz,idx=0,rnz,i,j,k,m;
251:   PetscInt       *row,*col;

254:   MatGetBlockSize(A,&bs);
255:   M = A->rmap->N/bs;
256:   *v = aa->a;
257:   if (reuse == MAT_INITIAL_MATRIX) {
258:     ai   = aa->i; aj = aa->j;
259:     nz   = bs2*aa->nz;
260:     *nnz = nz;
261:     PetscMalloc1(2*nz, &row);
262:     col  = row + nz;

264:     for (i=0; i<M; i++) {
265:       ajj = aj + ai[i];
266:       rnz = ai[i+1] - ai[i];
267:       for (k=0; k<rnz; k++) {
268:         for (j=0; j<bs; j++) {
269:           for (m=0; m<bs; m++) {
270:             row[idx]   = i*bs + m + shift;
271:             col[idx++] = bs*(ajj[k]) + j + shift;
272:           }
273:         }
274:       }
275:     }
276:     *r = row; *c = col;
277:   }
278:   return(0);
279: }

281: PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
282: {
283:   const PetscInt *ai, *aj,*ajj,M=A->rmap->n;
284:   PetscInt       nz,rnz,i,j;
286:   PetscInt       *row,*col;
287:   Mat_SeqSBAIJ   *aa=(Mat_SeqSBAIJ*)A->data;

290:   *v = aa->a;
291:   if (reuse == MAT_INITIAL_MATRIX) {
292:     nz   = aa->nz;
293:     ai   = aa->i;
294:     aj   = aa->j;
295:     *v   = aa->a;
296:     *nnz = nz;
297:     PetscMalloc1(2*nz, &row);
298:     col  = row + nz;

300:     nz = 0;
301:     for (i=0; i<M; i++) {
302:       rnz = ai[i+1] - ai[i];
303:       ajj = aj + ai[i];
304:       for (j=0; j<rnz; j++) {
305:         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
306:       }
307:     }
308:     *r = row; *c = col;
309:   }
310:   return(0);
311: }

313: PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
314: {
315:   const PetscInt    *ai,*aj,*ajj,*adiag,M=A->rmap->n;
316:   PetscInt          nz,rnz,i,j;
317:   const PetscScalar *av,*v1;
318:   PetscScalar       *val;
319:   PetscErrorCode    ierr;
320:   PetscInt          *row,*col;
321:   Mat_SeqAIJ        *aa=(Mat_SeqAIJ*)A->data;
322:   PetscBool         missing;

325:   ai   =aa->i; aj=aa->j;av=aa->a;
326:   adiag=aa->diag;
327:   MatMissingDiagonal_SeqAIJ(A,&missing,&i);
328:   if (reuse == MAT_INITIAL_MATRIX) {
329:     /* count nz in the uppper triangular part of A */
330:     nz = 0;
331:     if (missing) {
332:       for (i=0; i<M; i++) {
333:         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
334:           for (j=ai[i];j<ai[i+1];j++) {
335:             if (aj[j] < i) continue;
336:             nz++;
337:           }
338:         } else {
339:           nz += ai[i+1] - adiag[i];
340:         }
341:       }
342:     } else {
343:       for (i=0; i<M; i++) nz += ai[i+1] - adiag[i];
344:     }
345:     *nnz = nz;

347:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
348:     col  = row + nz;
349:     val  = (PetscScalar*)(col + nz);

351:     nz = 0;
352:     if (missing) {
353:       for (i=0; i<M; i++) {
354:         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
355:           for (j=ai[i];j<ai[i+1];j++) {
356:             if (aj[j] < i) continue;
357:             row[nz] = i+shift;
358:             col[nz] = aj[j]+shift;
359:             val[nz] = av[j];
360:             nz++;
361:           }
362:         } else {
363:           rnz = ai[i+1] - adiag[i];
364:           ajj = aj + adiag[i];
365:           v1  = av + adiag[i];
366:           for (j=0; j<rnz; j++) {
367:             row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
368:           }
369:         }
370:       }
371:     } else {
372:       for (i=0; i<M; i++) {
373:         rnz = ai[i+1] - adiag[i];
374:         ajj = aj + adiag[i];
375:         v1  = av + adiag[i];
376:         for (j=0; j<rnz; j++) {
377:           row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
378:         }
379:       }
380:     }
381:     *r = row; *c = col; *v = val;
382:   } else {
383:     nz = 0; val = *v;
384:     if (missing) {
385:       for (i=0; i <M; i++) {
386:         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
387:           for (j=ai[i];j<ai[i+1];j++) {
388:             if (aj[j] < i) continue;
389:             val[nz++] = av[j];
390:           }
391:         } else {
392:           rnz = ai[i+1] - adiag[i];
393:           v1  = av + adiag[i];
394:           for (j=0; j<rnz; j++) {
395:             val[nz++] = v1[j];
396:           }
397:         }
398:       }
399:     } else {
400:       for (i=0; i <M; i++) {
401:         rnz = ai[i+1] - adiag[i];
402:         v1  = av + adiag[i];
403:         for (j=0; j<rnz; j++) {
404:           val[nz++] = v1[j];
405:         }
406:       }
407:     }
408:   }
409:   return(0);
410: }

412: PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
413: {
414:   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
415:   PetscErrorCode    ierr;
416:   PetscInt          rstart,nz,i,j,jj,irow,countA,countB;
417:   PetscInt          *row,*col;
418:   const PetscScalar *av, *bv,*v1,*v2;
419:   PetscScalar       *val;
420:   Mat_MPISBAIJ      *mat = (Mat_MPISBAIJ*)A->data;
421:   Mat_SeqSBAIJ      *aa  = (Mat_SeqSBAIJ*)(mat->A)->data;
422:   Mat_SeqBAIJ       *bb  = (Mat_SeqBAIJ*)(mat->B)->data;

425:   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
426:   av=aa->a; bv=bb->a;

428:   garray = mat->garray;

430:   if (reuse == MAT_INITIAL_MATRIX) {
431:     nz   = aa->nz + bb->nz;
432:     *nnz = nz;
433:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
434:     col  = row + nz;
435:     val  = (PetscScalar*)(col + nz);

437:     *r = row; *c = col; *v = val;
438:   } else {
439:     row = *r; col = *c; val = *v;
440:   }

442:   jj = 0; irow = rstart;
443:   for (i=0; i<m; i++) {
444:     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
445:     countA = ai[i+1] - ai[i];
446:     countB = bi[i+1] - bi[i];
447:     bjj    = bj + bi[i];
448:     v1     = av + ai[i];
449:     v2     = bv + bi[i];

451:     /* A-part */
452:     for (j=0; j<countA; j++) {
453:       if (reuse == MAT_INITIAL_MATRIX) {
454:         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
455:       }
456:       val[jj++] = v1[j];
457:     }

459:     /* B-part */
460:     for (j=0; j < countB; j++) {
461:       if (reuse == MAT_INITIAL_MATRIX) {
462:         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
463:       }
464:       val[jj++] = v2[j];
465:     }
466:     irow++;
467:   }
468:   return(0);
469: }

471: PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
472: {
473:   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
474:   PetscErrorCode    ierr;
475:   PetscInt          rstart,nz,i,j,jj,irow,countA,countB;
476:   PetscInt          *row,*col;
477:   const PetscScalar *av, *bv,*v1,*v2;
478:   PetscScalar       *val;
479:   Mat_MPIAIJ        *mat = (Mat_MPIAIJ*)A->data;
480:   Mat_SeqAIJ        *aa  = (Mat_SeqAIJ*)(mat->A)->data;
481:   Mat_SeqAIJ        *bb  = (Mat_SeqAIJ*)(mat->B)->data;

484:   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
485:   av=aa->a; bv=bb->a;

487:   garray = mat->garray;

489:   if (reuse == MAT_INITIAL_MATRIX) {
490:     nz   = aa->nz + bb->nz;
491:     *nnz = nz;
492:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
493:     col  = row + nz;
494:     val  = (PetscScalar*)(col + nz);

496:     *r = row; *c = col; *v = val;
497:   } else {
498:     row = *r; col = *c; val = *v;
499:   }

501:   jj = 0; irow = rstart;
502:   for (i=0; i<m; i++) {
503:     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
504:     countA = ai[i+1] - ai[i];
505:     countB = bi[i+1] - bi[i];
506:     bjj    = bj + bi[i];
507:     v1     = av + ai[i];
508:     v2     = bv + bi[i];

510:     /* A-part */
511:     for (j=0; j<countA; j++) {
512:       if (reuse == MAT_INITIAL_MATRIX) {
513:         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
514:       }
515:       val[jj++] = v1[j];
516:     }

518:     /* B-part */
519:     for (j=0; j < countB; j++) {
520:       if (reuse == MAT_INITIAL_MATRIX) {
521:         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
522:       }
523:       val[jj++] = v2[j];
524:     }
525:     irow++;
526:   }
527:   return(0);
528: }

530: PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
531: {
532:   Mat_MPIBAIJ       *mat    = (Mat_MPIBAIJ*)A->data;
533:   Mat_SeqBAIJ       *aa     = (Mat_SeqBAIJ*)(mat->A)->data;
534:   Mat_SeqBAIJ       *bb     = (Mat_SeqBAIJ*)(mat->B)->data;
535:   const PetscInt    *ai     = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj;
536:   const PetscInt    *garray = mat->garray,mbs=mat->mbs,rstart=A->rmap->rstart;
537:   const PetscInt    bs2=mat->bs2;
538:   PetscErrorCode    ierr;
539:   PetscInt          bs,nz,i,j,k,n,jj,irow,countA,countB,idx;
540:   PetscInt          *row,*col;
541:   const PetscScalar *av=aa->a, *bv=bb->a,*v1,*v2;
542:   PetscScalar       *val;

545:   MatGetBlockSize(A,&bs);
546:   if (reuse == MAT_INITIAL_MATRIX) {
547:     nz   = bs2*(aa->nz + bb->nz);
548:     *nnz = nz;
549:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
550:     col  = row + nz;
551:     val  = (PetscScalar*)(col + nz);

553:     *r = row; *c = col; *v = val;
554:   } else {
555:     row = *r; col = *c; val = *v;
556:   }

558:   jj = 0; irow = rstart;
559:   for (i=0; i<mbs; i++) {
560:     countA = ai[i+1] - ai[i];
561:     countB = bi[i+1] - bi[i];
562:     ajj    = aj + ai[i];
563:     bjj    = bj + bi[i];
564:     v1     = av + bs2*ai[i];
565:     v2     = bv + bs2*bi[i];

567:     idx = 0;
568:     /* A-part */
569:     for (k=0; k<countA; k++) {
570:       for (j=0; j<bs; j++) {
571:         for (n=0; n<bs; n++) {
572:           if (reuse == MAT_INITIAL_MATRIX) {
573:             row[jj] = irow + n + shift;
574:             col[jj] = rstart + bs*ajj[k] + j + shift;
575:           }
576:           val[jj++] = v1[idx++];
577:         }
578:       }
579:     }

581:     idx = 0;
582:     /* B-part */
583:     for (k=0; k<countB; k++) {
584:       for (j=0; j<bs; j++) {
585:         for (n=0; n<bs; n++) {
586:           if (reuse == MAT_INITIAL_MATRIX) {
587:             row[jj] = irow + n + shift;
588:             col[jj] = bs*garray[bjj[k]] + j + shift;
589:           }
590:           val[jj++] = v2[idx++];
591:         }
592:       }
593:     }
594:     irow += bs;
595:   }
596:   return(0);
597: }

599: PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
600: {
601:   const PetscInt    *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
602:   PetscErrorCode    ierr;
603:   PetscInt          rstart,nz,nza,nzb,i,j,jj,irow,countA,countB;
604:   PetscInt          *row,*col;
605:   const PetscScalar *av, *bv,*v1,*v2;
606:   PetscScalar       *val;
607:   Mat_MPIAIJ        *mat =  (Mat_MPIAIJ*)A->data;
608:   Mat_SeqAIJ        *aa  =(Mat_SeqAIJ*)(mat->A)->data;
609:   Mat_SeqAIJ        *bb  =(Mat_SeqAIJ*)(mat->B)->data;

612:   ai=aa->i; aj=aa->j; adiag=aa->diag;
613:   bi=bb->i; bj=bb->j; garray = mat->garray;
614:   av=aa->a; bv=bb->a;

616:   rstart = A->rmap->rstart;

618:   if (reuse == MAT_INITIAL_MATRIX) {
619:     nza = 0;    /* num of upper triangular entries in mat->A, including diagonals */
620:     nzb = 0;    /* num of upper triangular entries in mat->B */
621:     for (i=0; i<m; i++) {
622:       nza   += (ai[i+1] - adiag[i]);
623:       countB = bi[i+1] - bi[i];
624:       bjj    = bj + bi[i];
625:       for (j=0; j<countB; j++) {
626:         if (garray[bjj[j]] > rstart) nzb++;
627:       }
628:     }

630:     nz   = nza + nzb; /* total nz of upper triangular part of mat */
631:     *nnz = nz;
632:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
633:     col  = row + nz;
634:     val  = (PetscScalar*)(col + nz);

636:     *r = row; *c = col; *v = val;
637:   } else {
638:     row = *r; col = *c; val = *v;
639:   }

641:   jj = 0; irow = rstart;
642:   for (i=0; i<m; i++) {
643:     ajj    = aj + adiag[i];                 /* ptr to the beginning of the diagonal of this row */
644:     v1     = av + adiag[i];
645:     countA = ai[i+1] - adiag[i];
646:     countB = bi[i+1] - bi[i];
647:     bjj    = bj + bi[i];
648:     v2     = bv + bi[i];

650:     /* A-part */
651:     for (j=0; j<countA; j++) {
652:       if (reuse == MAT_INITIAL_MATRIX) {
653:         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
654:       }
655:       val[jj++] = v1[j];
656:     }

658:     /* B-part */
659:     for (j=0; j < countB; j++) {
660:       if (garray[bjj[j]] > rstart) {
661:         if (reuse == MAT_INITIAL_MATRIX) {
662:           row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
663:         }
664:         val[jj++] = v2[j];
665:       }
666:     }
667:     irow++;
668:   }
669:   return(0);
670: }

672: PetscErrorCode MatDestroy_MUMPS(Mat A)
673: {
674:   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->data;

678:   PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);
679:   VecScatterDestroy(&mumps->scat_rhs);
680:   VecScatterDestroy(&mumps->scat_sol);
681:   VecDestroy(&mumps->b_seq);
682:   VecDestroy(&mumps->x_seq);
683:   PetscFree(mumps->id.perm_in);
684:   PetscFree(mumps->irn);
685:   PetscFree(mumps->info);
686:   MatMumpsResetSchur_Private(mumps);
687:   mumps->id.job = JOB_END;
688:   PetscMUMPS_c(&mumps->id);
689:   MPI_Comm_free(&mumps->comm_mumps);
690:   PetscFree(A->data);

692:   /* clear composed functions */
693:   PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverType_C",NULL);
694:   PetscObjectComposeFunction((PetscObject)A,"MatFactorSetSchurIS_C",NULL);
695:   PetscObjectComposeFunction((PetscObject)A,"MatFactorCreateSchurComplement_C",NULL);
696:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetIcntl_C",NULL);
697:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetIcntl_C",NULL);
698:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetCntl_C",NULL);
699:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetCntl_C",NULL);
700:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfo_C",NULL);
701:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfog_C",NULL);
702:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfo_C",NULL);
703:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfog_C",NULL);
704:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInverse_C",NULL);
705:   return(0);
706: }

708: PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x)
709: {
710:   Mat_MUMPS        *mumps=(Mat_MUMPS*)A->data;
711:   PetscScalar      *array;
712:   Vec              b_seq;
713:   IS               is_iden,is_petsc;
714:   PetscErrorCode   ierr;
715:   PetscInt         i;
716:   PetscBool        second_solve = PETSC_FALSE;
717:   static PetscBool cite1 = PETSC_FALSE,cite2 = PETSC_FALSE;

720:   PetscCitationsRegister("@article{MUMPS01,\n  author = {P.~R. Amestoy and I.~S. Duff and J.-Y. L'Excellent and J. Koster},\n  title = {A fully asynchronous multifrontal solver using distributed dynamic scheduling},\n  journal = {SIAM Journal on Matrix Analysis and Applications},\n  volume = {23},\n  number = {1},\n  pages = {15--41},\n  year = {2001}\n}\n",&cite1);
721:   PetscCitationsRegister("@article{MUMPS02,\n  author = {P.~R. Amestoy and A. Guermouche and J.-Y. L'Excellent and S. Pralet},\n  title = {Hybrid scheduling for the parallel solution of linear systems},\n  journal = {Parallel Computing},\n  volume = {32},\n  number = {2},\n  pages = {136--156},\n  year = {2006}\n}\n",&cite2);

723:   if (A->factorerrortype) {
724:     PetscInfo2(A,"MatSolve is called with singular matrix factor, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
725:     VecSetInf(x);
726:     return(0);
727:   }

729:   mumps->id.nrhs = 1;
730:   b_seq          = mumps->b_seq;
731:   if (mumps->size > 1) {
732:     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
733:     VecScatterBegin(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);
734:     VecScatterEnd(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);
735:     if (!mumps->myid) {VecGetArray(b_seq,&array);}
736:   } else {  /* size == 1 */
737:     VecCopy(b,x);
738:     VecGetArray(x,&array);
739:   }
740:   if (!mumps->myid) { /* define rhs on the host */
741:     mumps->id.nrhs = 1;
742:     mumps->id.rhs = (MumpsScalar*)array;
743:   }

745:   /*
746:      handle condensation step of Schur complement (if any)
747:      We set by default ICNTL(26) == -1 when Schur indices have been provided by the user.
748:      According to MUMPS (5.0.0) manual, any value should be harmful during the factorization phase
749:      Unless the user provides a valid value for ICNTL(26), MatSolve and MatMatSolve routines solve the full system.
750:      This requires an extra call to PetscMUMPS_c and the computation of the factors for S
751:   */
752:   if (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2) {
753:     if (mumps->size > 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Parallel Schur complements not yet supported from PETSc\n");
754:     second_solve = PETSC_TRUE;
755:     MatMumpsHandleSchur_Private(A,PETSC_FALSE);
756:   }
757:   /* solve phase */
758:   /*-------------*/
759:   mumps->id.job = JOB_SOLVE;
760:   PetscMUMPS_c(&mumps->id);
761:   if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",mumps->id.INFOG(1));

763:   /* handle expansion step of Schur complement (if any) */
764:   if (second_solve) {
765:     MatMumpsHandleSchur_Private(A,PETSC_TRUE);
766:   }

768:   if (mumps->size > 1) { /* convert mumps distributed solution to petsc mpi x */
769:     if (mumps->scat_sol && mumps->ICNTL9_pre != mumps->id.ICNTL(9)) {
770:       /* when id.ICNTL(9) changes, the contents of lsol_loc may change (not its size, lsol_loc), recreates scat_sol */
771:       VecScatterDestroy(&mumps->scat_sol);
772:     }
773:     if (!mumps->scat_sol) { /* create scatter scat_sol */
774:       ISCreateStride(PETSC_COMM_SELF,mumps->id.lsol_loc,0,1,&is_iden); /* from */
775:       for (i=0; i<mumps->id.lsol_loc; i++) {
776:         mumps->id.isol_loc[i] -= 1; /* change Fortran style to C style */
777:       }
778:       ISCreateGeneral(PETSC_COMM_SELF,mumps->id.lsol_loc,mumps->id.isol_loc,PETSC_COPY_VALUES,&is_petsc);  /* to */
779:       VecScatterCreate(mumps->x_seq,is_iden,x,is_petsc,&mumps->scat_sol);
780:       ISDestroy(&is_iden);
781:       ISDestroy(&is_petsc);

783:       mumps->ICNTL9_pre = mumps->id.ICNTL(9); /* save current value of id.ICNTL(9) */
784:     }

786:     VecScatterBegin(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);
787:     VecScatterEnd(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);
788:   }
789:   return(0);
790: }

792: PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x)
793: {
794:   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->data;

798:   mumps->id.ICNTL(9) = 0;
799:   MatSolve_MUMPS(A,b,x);
800:   mumps->id.ICNTL(9) = 1;
801:   return(0);
802: }

804: PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X)
805: {
807:   Mat            Bt = NULL;
808:   PetscBool      flg, flgT;
809:   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->data;
810:   PetscInt       i,nrhs,M;
811:   PetscScalar    *array,*bray;

814:   PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
815:   PetscObjectTypeCompare((PetscObject)B,MATTRANSPOSEMAT,&flgT);
816:   if (flgT) {
817:     if (mumps->size > 1) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
818:     MatTransposeGetMat(B,&Bt);
819:   } else {
820:     if (!flg) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
821:   }

823:   MatGetSize(B,&M,&nrhs);
824:   mumps->id.nrhs = nrhs;
825:   mumps->id.lrhs = M;

827:   PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
828:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");

830:   if (B->rmap->n != X->rmap->n) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B and X must have same row distribution");

832:   if (mumps->size == 1) {
833:     PetscScalar *aa;
834:     PetscInt    spnr,*ia,*ja;
835:     PetscBool   second_solve = PETSC_FALSE;

837:     /* copy B to X */
838:     MatDenseGetArray(X,&array);
839:     mumps->id.rhs = (MumpsScalar*)array;
840:     if (!Bt) {
841:       MatDenseGetArray(B,&bray);
842:       PetscMemcpy(array,bray,M*nrhs*sizeof(PetscScalar));
843:       MatDenseRestoreArray(B,&bray);
844:     } else {
845:       PetscBool done;

847:       MatSeqAIJGetArray(Bt,&aa);
848:       MatGetRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&done);
849:       if (!done) SETERRQ(PetscObjectComm((PetscObject)Bt),PETSC_ERR_ARG_WRONG,"Cannot get IJ structure");
850:       mumps->id.irhs_ptr    = ia;
851:       mumps->id.irhs_sparse = ja;
852:       mumps->id.nz_rhs      = ia[spnr] - 1;
853:       mumps->id.rhs_sparse  = (MumpsScalar*)aa;
854:       mumps->id.ICNTL(20)   = 1;
855:     }
856:     /* handle condensation step of Schur complement (if any) */
857:     if (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2) {
858:       second_solve = PETSC_TRUE;
859:       MatMumpsHandleSchur_Private(A,PETSC_FALSE);
860:     }
861:     /* solve phase */
862:     /*-------------*/
863:     mumps->id.job = JOB_SOLVE;
864:     PetscMUMPS_c(&mumps->id);
865:     if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",mumps->id.INFOG(1));

867:     /* handle expansion step of Schur complement (if any) */
868:     if (second_solve) {
869:       MatMumpsHandleSchur_Private(A,PETSC_TRUE);
870:     }
871:     if (Bt) {
872:       PetscBool done;

874:       MatSeqAIJRestoreArray(Bt,&aa);
875:       MatRestoreRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&done);
876:       if (!done) SETERRQ(PetscObjectComm((PetscObject)Bt),PETSC_ERR_ARG_WRONG,"Cannot restore IJ structure");
877:       mumps->id.ICNTL(20) = 0;
878:     }
879:     MatDenseRestoreArray(X,&array);
880:   } else {  /*--------- parallel case --------*/
881:     PetscInt       lsol_loc,nlsol_loc,*isol_loc,*idx,*iidx,*idxx,*isol_loc_save;
882:     MumpsScalar    *sol_loc,*sol_loc_save;
883:     IS             is_to,is_from;
884:     PetscInt       k,proc,j,m;
885:     const PetscInt *rstart;
886:     Vec            v_mpi,b_seq,x_seq;
887:     VecScatter     scat_rhs,scat_sol;

889:     if (mumps->size > 1 && mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Parallel Schur complements not yet supported from PETSc\n");

891:     /* create x_seq to hold local solution */
892:     isol_loc_save = mumps->id.isol_loc; /* save it for MatSovle() */
893:     sol_loc_save  = mumps->id.sol_loc;

895:     lsol_loc  = mumps->id.INFO(23);
896:     nlsol_loc = nrhs*lsol_loc;     /* length of sol_loc */
897:     PetscMalloc2(nlsol_loc,&sol_loc,nlsol_loc,&isol_loc);
898:     mumps->id.sol_loc = (MumpsScalar*)sol_loc;
899:     mumps->id.isol_loc = isol_loc;

901:     VecCreateSeqWithArray(PETSC_COMM_SELF,1,nlsol_loc,(PetscScalar*)sol_loc,&x_seq);

903:     /* copy rhs matrix B into vector v_mpi */
904:     MatGetLocalSize(B,&m,NULL);
905:     MatDenseGetArray(B,&bray);
906:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)B),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi);
907:     MatDenseRestoreArray(B,&bray);

909:     /* scatter v_mpi to b_seq because MUMPS only supports centralized rhs */
910:     /* idx: maps from k-th index of v_mpi to (i,j)-th global entry of B;
911:       iidx: inverse of idx, will be used by scattering xx_seq -> X       */
912:     PetscMalloc2(nrhs*M,&idx,nrhs*M,&iidx);
913:     MatGetOwnershipRanges(B,&rstart);
914:     k = 0;
915:     for (proc=0; proc<mumps->size; proc++){
916:       for (j=0; j<nrhs; j++){
917:         for (i=rstart[proc]; i<rstart[proc+1]; i++){
918:           iidx[j*M + i] = k;
919:           idx[k++]      = j*M + i;
920:         }
921:       }
922:     }

924:     if (!mumps->myid) {
925:       VecCreateSeq(PETSC_COMM_SELF,nrhs*M,&b_seq);
926:       ISCreateGeneral(PETSC_COMM_SELF,nrhs*M,idx,PETSC_COPY_VALUES,&is_to);
927:       ISCreateStride(PETSC_COMM_SELF,nrhs*M,0,1,&is_from);
928:     } else {
929:       VecCreateSeq(PETSC_COMM_SELF,0,&b_seq);
930:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_to);
931:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_from);
932:     }
933:     VecScatterCreate(v_mpi,is_from,b_seq,is_to,&scat_rhs);
934:     VecScatterBegin(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);
935:     ISDestroy(&is_to);
936:     ISDestroy(&is_from);
937:     VecScatterEnd(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);

939:     if (!mumps->myid) { /* define rhs on the host */
940:       VecGetArray(b_seq,&bray);
941:       mumps->id.rhs = (MumpsScalar*)bray;
942:       VecRestoreArray(b_seq,&bray);
943:     }

945:     /* solve phase */
946:     /*-------------*/
947:     mumps->id.job = JOB_SOLVE;
948:     PetscMUMPS_c(&mumps->id);
949:     if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",mumps->id.INFOG(1));

951:     /* scatter mumps distributed solution to petsc vector v_mpi, which shares local arrays with solution matrix X */
952:     MatDenseGetArray(X,&array);
953:     VecPlaceArray(v_mpi,array);

955:     /* create scatter scat_sol */
956:     PetscMalloc1(nlsol_loc,&idxx);
957:     ISCreateStride(PETSC_COMM_SELF,nlsol_loc,0,1,&is_from);
958:     for (i=0; i<lsol_loc; i++) {
959:       isol_loc[i] -= 1; /* change Fortran style to C style */
960:       idxx[i] = iidx[isol_loc[i]];
961:       for (j=1; j<nrhs; j++){
962:         idxx[j*lsol_loc+i] = iidx[isol_loc[i]+j*M];
963:       }
964:     }
965:     ISCreateGeneral(PETSC_COMM_SELF,nlsol_loc,idxx,PETSC_COPY_VALUES,&is_to);
966:     VecScatterCreate(x_seq,is_from,v_mpi,is_to,&scat_sol);
967:     VecScatterBegin(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);
968:     ISDestroy(&is_from);
969:     ISDestroy(&is_to);
970:     VecScatterEnd(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);
971:     MatDenseRestoreArray(X,&array);

973:     /* free spaces */
974:     mumps->id.sol_loc = sol_loc_save;
975:     mumps->id.isol_loc = isol_loc_save;

977:     PetscFree2(sol_loc,isol_loc);
978:     PetscFree2(idx,iidx);
979:     PetscFree(idxx);
980:     VecDestroy(&x_seq);
981:     VecDestroy(&v_mpi);
982:     VecDestroy(&b_seq);
983:     VecScatterDestroy(&scat_rhs);
984:     VecScatterDestroy(&scat_sol);
985:   }
986:   return(0);
987: }

989: #if !defined(PETSC_USE_COMPLEX)
990: /*
991:   input:
992:    F:        numeric factor
993:   output:
994:    nneg:     total number of negative pivots
995:    nzero:    total number of zero pivots
996:    npos:     (global dimension of F) - nneg - nzero
997: */
998: PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
999: {
1000:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1002:   PetscMPIInt    size;

1005:   MPI_Comm_size(PetscObjectComm((PetscObject)F),&size);
1006:   /* 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 */
1007:   if (size > 1 && mumps->id.ICNTL(13) != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"ICNTL(13)=%d. -mat_mumps_icntl_13 must be set as 1 for correct global matrix inertia\n",mumps->id.INFOG(13));

1009:   if (nneg) *nneg = mumps->id.INFOG(12);
1010:   if (nzero || npos) {
1011:     if (mumps->id.ICNTL(24) != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"-mat_mumps_icntl_24 must be set as 1 for null pivot row detection");
1012:     if (nzero) *nzero = mumps->id.INFOG(28);
1013:     if (npos) *npos   = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28));
1014:   }
1015:   return(0);
1016: }
1017: #endif

1019: PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
1020: {
1021:   Mat_MUMPS      *mumps =(Mat_MUMPS*)(F)->data;
1023:   PetscBool      isMPIAIJ;

1026:   if (mumps->id.INFOG(1) < 0) {
1027:     if (mumps->id.INFOG(1) == -6) {
1028:       PetscInfo2(A,"MatFactorNumeric is called with singular matrix structure, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1029:     }
1030:     PetscInfo2(A,"MatFactorNumeric is called after analysis phase fails, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1031:     return(0);
1032:   }

1034:   (*mumps->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);

1036:   /* numerical factorization phase */
1037:   /*-------------------------------*/
1038:   mumps->id.job = JOB_FACTNUMERIC;
1039:   if (!mumps->id.ICNTL(18)) { /* A is centralized */
1040:     if (!mumps->myid) {
1041:       mumps->id.a = (MumpsScalar*)mumps->val;
1042:     }
1043:   } else {
1044:     mumps->id.a_loc = (MumpsScalar*)mumps->val;
1045:   }
1046:   PetscMUMPS_c(&mumps->id);
1047:   if (mumps->id.INFOG(1) < 0) {
1048:     if (A->erroriffailure) {
1049:       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1050:     } else {
1051:       if (mumps->id.INFOG(1) == -10) { /* numerically singular matrix */
1052:         PetscInfo2(F,"matrix is numerically singular, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1053:         F->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1054:       } else if (mumps->id.INFOG(1) == -13) {
1055:         PetscInfo2(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, cannot allocate required memory %d megabytes\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1056:         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1057:       } else if (mumps->id.INFOG(1) == -8 || mumps->id.INFOG(1) == -9 || (-16 < mumps->id.INFOG(1) && mumps->id.INFOG(1) < -10) ) {
1058:         PetscInfo2(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d, problem with workarray \n",mumps->id.INFOG(1),mumps->id.INFO(2));
1059:         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1060:       } else {
1061:         PetscInfo2(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1062:         F->factorerrortype = MAT_FACTOR_OTHER;
1063:       }
1064:     }
1065:   }
1066:   if (!mumps->myid && mumps->id.ICNTL(16) > 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"  mumps->id.ICNTL(16):=%d\n",mumps->id.INFOG(16));

1068:   F->assembled    = PETSC_TRUE;
1069:   mumps->matstruc = SAME_NONZERO_PATTERN;
1070:   if (F->schur) { /* reset Schur status to unfactored */
1071:     if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
1072:       mumps->id.ICNTL(19) = 2;
1073:       MatTranspose(F->schur,MAT_INPLACE_MATRIX,&F->schur);
1074:     }
1075:     MatFactorRestoreSchurComplement(F,NULL,MAT_FACTOR_SCHUR_UNFACTORED);
1076:   }

1078:   /* just to be sure that ICNTL(19) value returned by a call from MatMumpsGetIcntl is always consistent */
1079:   if (!mumps->sym && mumps->id.ICNTL(19) && mumps->id.ICNTL(19) != 1) mumps->id.ICNTL(19) = 3;

1081:   if (mumps->size > 1) {
1082:     PetscInt    lsol_loc;
1083:     PetscScalar *sol_loc;

1085:     PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);

1087:     /* distributed solution; Create x_seq=sol_loc for repeated use */
1088:     if (mumps->x_seq) {
1089:       VecScatterDestroy(&mumps->scat_sol);
1090:       PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);
1091:       VecDestroy(&mumps->x_seq);
1092:     }
1093:     lsol_loc = mumps->id.INFO(23); /* length of sol_loc */
1094:     PetscMalloc2(lsol_loc,&sol_loc,lsol_loc,&mumps->id.isol_loc);
1095:     mumps->id.lsol_loc = lsol_loc;
1096:     mumps->id.sol_loc = (MumpsScalar*)sol_loc;
1097:     VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq);
1098:   }
1099:   return(0);
1100: }

1102: /* Sets MUMPS options from the options database */
1103: PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A)
1104: {
1105:   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1107:   PetscInt       icntl,info[40],i,ninfo=40;
1108:   PetscBool      flg;

1111:   PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MUMPS Options","Mat");
1112:   PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);
1113:   if (flg) mumps->id.ICNTL(1) = icntl;
1114:   PetscOptionsInt("-mat_mumps_icntl_2","ICNTL(2): output stream for diagnostic printing, statistics, and warning","None",mumps->id.ICNTL(2),&icntl,&flg);
1115:   if (flg) mumps->id.ICNTL(2) = icntl;
1116:   PetscOptionsInt("-mat_mumps_icntl_3","ICNTL(3): output stream for global information, collected on the host","None",mumps->id.ICNTL(3),&icntl,&flg);
1117:   if (flg) mumps->id.ICNTL(3) = icntl;

1119:   PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",mumps->id.ICNTL(4),&icntl,&flg);
1120:   if (flg) mumps->id.ICNTL(4) = icntl;
1121:   if (mumps->id.ICNTL(4) || PetscLogPrintInfo) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */

1123:   PetscOptionsInt("-mat_mumps_icntl_6","ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7)","None",mumps->id.ICNTL(6),&icntl,&flg);
1124:   if (flg) mumps->id.ICNTL(6) = icntl;

1126:   PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis","None",mumps->id.ICNTL(7),&icntl,&flg);
1127:   if (flg) {
1128:     if (icntl== 1 && mumps->size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n");
1129:     else mumps->id.ICNTL(7) = icntl;
1130:   }

1132:   PetscOptionsInt("-mat_mumps_icntl_8","ICNTL(8): scaling strategy (-2 to 8 or 77)","None",mumps->id.ICNTL(8),&mumps->id.ICNTL(8),NULL);
1133:   /* PetscOptionsInt("-mat_mumps_icntl_9","ICNTL(9): computes the solution using A or A^T","None",mumps->id.ICNTL(9),&mumps->id.ICNTL(9),NULL); handled by MatSolveTranspose_MUMPS() */
1134:   PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",mumps->id.ICNTL(10),&mumps->id.ICNTL(10),NULL);
1135:   PetscOptionsInt("-mat_mumps_icntl_11","ICNTL(11): statistics related to an error analysis (via -ksp_view)","None",mumps->id.ICNTL(11),&mumps->id.ICNTL(11),NULL);
1136:   PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): an ordering strategy for symmetric matrices (0 to 3)","None",mumps->id.ICNTL(12),&mumps->id.ICNTL(12),NULL);
1137:   PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting","None",mumps->id.ICNTL(13),&mumps->id.ICNTL(13),NULL);
1138:   PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage increase in the estimated working space","None",mumps->id.ICNTL(14),&mumps->id.ICNTL(14),NULL);
1139:   PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): computes the Schur complement","None",mumps->id.ICNTL(19),&mumps->id.ICNTL(19),NULL);
1140:   if (mumps->id.ICNTL(19) <= 0 || mumps->id.ICNTL(19) > 3) { /* reset any schur data (if any) */
1141:     MatDestroy(&F->schur);
1142:     MatMumpsResetSchur_Private(mumps);
1143:   }
1144:   /* PetscOptionsInt("-mat_mumps_icntl_20","ICNTL(20): the format (dense or sparse) of the right-hand sides","None",mumps->id.ICNTL(20),&mumps->id.ICNTL(20),NULL); -- sparse rhs is not supported in PETSc API */
1145:   /* PetscOptionsInt("-mat_mumps_icntl_21","ICNTL(21): the distribution (centralized or distributed) of the solution vectors","None",mumps->id.ICNTL(21),&mumps->id.ICNTL(21),NULL); we only use distributed solution vector */

1147:   PetscOptionsInt("-mat_mumps_icntl_22","ICNTL(22): in-core/out-of-core factorization and solve (0 or 1)","None",mumps->id.ICNTL(22),&mumps->id.ICNTL(22),NULL);
1148:   PetscOptionsInt("-mat_mumps_icntl_23","ICNTL(23): max size of the working memory (MB) that can allocate per processor","None",mumps->id.ICNTL(23),&mumps->id.ICNTL(23),NULL);
1149:   PetscOptionsInt("-mat_mumps_icntl_24","ICNTL(24): detection of null pivot rows (0 or 1)","None",mumps->id.ICNTL(24),&mumps->id.ICNTL(24),NULL);
1150:   if (mumps->id.ICNTL(24)) {
1151:     mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */
1152:   }

1154:   PetscOptionsInt("-mat_mumps_icntl_25","ICNTL(25): computes a solution of a deficient matrix and a null space basis","None",mumps->id.ICNTL(25),&mumps->id.ICNTL(25),NULL);
1155:   PetscOptionsInt("-mat_mumps_icntl_26","ICNTL(26): drives the solution phase if a Schur complement matrix","None",mumps->id.ICNTL(26),&mumps->id.ICNTL(26),NULL);
1156:   PetscOptionsInt("-mat_mumps_icntl_27","ICNTL(27): the blocking size for multiple right-hand sides","None",mumps->id.ICNTL(27),&mumps->id.ICNTL(27),NULL);
1157:   PetscOptionsInt("-mat_mumps_icntl_28","ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering","None",mumps->id.ICNTL(28),&mumps->id.ICNTL(28),NULL);
1158:   PetscOptionsInt("-mat_mumps_icntl_29","ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis","None",mumps->id.ICNTL(29),&mumps->id.ICNTL(29),NULL);
1159:   /* PetscOptionsInt("-mat_mumps_icntl_30","ICNTL(30): compute user-specified set of entries in inv(A)","None",mumps->id.ICNTL(30),&mumps->id.ICNTL(30),NULL); */ /* call MatMumpsGetInverse() directly */
1160:   PetscOptionsInt("-mat_mumps_icntl_31","ICNTL(31): indicates which factors may be discarded during factorization","None",mumps->id.ICNTL(31),&mumps->id.ICNTL(31),NULL);
1161:   /* PetscOptionsInt("-mat_mumps_icntl_32","ICNTL(32): performs the forward elemination of the right-hand sides during factorization","None",mumps->id.ICNTL(32),&mumps->id.ICNTL(32),NULL);  -- not supported by PETSc API */
1162:   PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),NULL);
1163:   PetscOptionsInt("-mat_mumps_icntl_35","ICNTL(35): activates Block Lock Rank (BLR) based factorization","None",mumps->id.ICNTL(35),&mumps->id.ICNTL(35),NULL);

1165:   PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),NULL);
1166:   PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",mumps->id.CNTL(2),&mumps->id.CNTL(2),NULL);
1167:   PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),NULL);
1168:   PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",mumps->id.CNTL(4),&mumps->id.CNTL(4),NULL);
1169:   PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",mumps->id.CNTL(5),&mumps->id.CNTL(5),NULL);
1170:   PetscOptionsReal("-mat_mumps_cntl_7","CNTL(7): dropping parameter used during BLR","None",mumps->id.CNTL(7),&mumps->id.CNTL(7),NULL);

1172:   PetscOptionsString("-mat_mumps_ooc_tmpdir", "out of core directory", "None", mumps->id.ooc_tmpdir, mumps->id.ooc_tmpdir, 256, NULL);

1174:   PetscOptionsIntArray("-mat_mumps_view_info","request INFO local to each processor","",info,&ninfo,NULL);
1175:   if (ninfo) {
1176:     if (ninfo > 40) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"number of INFO %d must <= 40\n",ninfo);
1177:     PetscMalloc1(ninfo,&mumps->info);
1178:     mumps->ninfo = ninfo;
1179:     for (i=0; i<ninfo; i++) {
1180:       if (info[i] < 0 || info[i]>40) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index of INFO %d must between 1 and 40\n",ninfo);
1181:       else  mumps->info[i] = info[i];
1182:     }
1183:   }

1185:   PetscOptionsEnd();
1186:   return(0);
1187: }

1189: PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps)
1190: {

1194:   MPI_Comm_rank(PetscObjectComm((PetscObject)A), &mumps->myid);
1195:   MPI_Comm_size(PetscObjectComm((PetscObject)A),&mumps->size);
1196:   MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(mumps->comm_mumps));

1198:   mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps);

1200:   mumps->id.job = JOB_INIT;
1201:   mumps->id.par = 1;  /* host participates factorizaton and solve */
1202:   mumps->id.sym = mumps->sym;
1203:   PetscMUMPS_c(&mumps->id);

1205:   mumps->scat_rhs     = NULL;
1206:   mumps->scat_sol     = NULL;

1208:   /* set PETSc-MUMPS default options - override MUMPS default */
1209:   mumps->id.ICNTL(3) = 0;
1210:   mumps->id.ICNTL(4) = 0;
1211:   if (mumps->size == 1) {
1212:     mumps->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
1213:   } else {
1214:     mumps->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
1215:     mumps->id.ICNTL(20) = 0;   /* rhs is in dense format */
1216:     mumps->id.ICNTL(21) = 1;   /* distributed solution */
1217:   }

1219:   /* schur */
1220:   mumps->id.size_schur      = 0;
1221:   mumps->id.listvar_schur   = NULL;
1222:   mumps->id.schur           = NULL;
1223:   mumps->sizeredrhs         = 0;
1224:   mumps->schur_sol          = NULL;
1225:   mumps->schur_sizesol      = 0;
1226:   return(0);
1227: }

1229: PetscErrorCode MatFactorSymbolic_MUMPS_ReportIfError(Mat F,Mat A,const MatFactorInfo *info,Mat_MUMPS *mumps)
1230: {

1234:   if (mumps->id.INFOG(1) < 0) {
1235:     if (A->erroriffailure) {
1236:       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
1237:     } else {
1238:       if (mumps->id.INFOG(1) == -6) {
1239:         PetscInfo2(F,"matrix is singular in structure, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1240:         F->factorerrortype = MAT_FACTOR_STRUCT_ZEROPIVOT;
1241:       } else if (mumps->id.INFOG(1) == -5 || mumps->id.INFOG(1) == -7) {
1242:         PetscInfo2(F,"problem of workspace, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1243:         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1244:       } else {
1245:         PetscInfo2(F,"Error reported by MUMPS in analysis phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1246:         F->factorerrortype = MAT_FACTOR_OTHER;
1247:       }
1248:     }
1249:   }
1250:   return(0);
1251: }

1253: /* Note Petsc r(=c) permutation is used when mumps->id.ICNTL(7)==1 with centralized assembled matrix input; otherwise r and c are ignored */
1254: PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1255: {
1256:   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1258:   Vec            b;
1259:   IS             is_iden;
1260:   const PetscInt M = A->rmap->N;

1263:   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;

1265:   /* Set MUMPS options from the options database */
1266:   PetscSetMUMPSFromOptions(F,A);

1268:   (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);

1270:   /* analysis phase */
1271:   /*----------------*/
1272:   mumps->id.job = JOB_FACTSYMBOLIC;
1273:   mumps->id.n   = M;
1274:   switch (mumps->id.ICNTL(18)) {
1275:   case 0:  /* centralized assembled matrix input */
1276:     if (!mumps->myid) {
1277:       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1278:       if (mumps->id.ICNTL(6)>1) {
1279:         mumps->id.a = (MumpsScalar*)mumps->val;
1280:       }
1281:       if (mumps->id.ICNTL(7) == 1) { /* use user-provide matrix ordering - assuming r = c ordering */
1282:         /*
1283:         PetscBool      flag;
1284:         ISEqual(r,c,&flag);
1285:         if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"row_perm != col_perm");
1286:         ISView(r,PETSC_VIEWER_STDOUT_SELF);
1287:          */
1288:         if (!mumps->myid) {
1289:           const PetscInt *idx;
1290:           PetscInt       i,*perm_in;

1292:           PetscMalloc1(M,&perm_in);
1293:           ISGetIndices(r,&idx);

1295:           mumps->id.perm_in = perm_in;
1296:           for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */
1297:           ISRestoreIndices(r,&idx);
1298:         }
1299:       }
1300:     }
1301:     break;
1302:   case 3:  /* distributed assembled matrix input (size>1) */
1303:     mumps->id.nz_loc = mumps->nz;
1304:     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1305:     if (mumps->id.ICNTL(6)>1) {
1306:       mumps->id.a_loc = (MumpsScalar*)mumps->val;
1307:     }
1308:     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1309:     if (!mumps->myid) {
1310:       VecCreateSeq(PETSC_COMM_SELF,A->rmap->N,&mumps->b_seq);
1311:       ISCreateStride(PETSC_COMM_SELF,A->rmap->N,0,1,&is_iden);
1312:     } else {
1313:       VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);
1314:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);
1315:     }
1316:     MatCreateVecs(A,NULL,&b);
1317:     VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);
1318:     ISDestroy(&is_iden);
1319:     VecDestroy(&b);
1320:     break;
1321:   }
1322:   PetscMUMPS_c(&mumps->id);
1323:   MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);

1325:   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1326:   F->ops->solve           = MatSolve_MUMPS;
1327:   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1328:   F->ops->matsolve        = MatMatSolve_MUMPS;
1329:   return(0);
1330: }

1332: /* Note the Petsc r and c permutations are ignored */
1333: PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1334: {
1335:   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1337:   Vec            b;
1338:   IS             is_iden;
1339:   const PetscInt M = A->rmap->N;

1342:   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;

1344:   /* Set MUMPS options from the options database */
1345:   PetscSetMUMPSFromOptions(F,A);

1347:   (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);

1349:   /* analysis phase */
1350:   /*----------------*/
1351:   mumps->id.job = JOB_FACTSYMBOLIC;
1352:   mumps->id.n   = M;
1353:   switch (mumps->id.ICNTL(18)) {
1354:   case 0:  /* centralized assembled matrix input */
1355:     if (!mumps->myid) {
1356:       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1357:       if (mumps->id.ICNTL(6)>1) {
1358:         mumps->id.a = (MumpsScalar*)mumps->val;
1359:       }
1360:     }
1361:     break;
1362:   case 3:  /* distributed assembled matrix input (size>1) */
1363:     mumps->id.nz_loc = mumps->nz;
1364:     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1365:     if (mumps->id.ICNTL(6)>1) {
1366:       mumps->id.a_loc = (MumpsScalar*)mumps->val;
1367:     }
1368:     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1369:     if (!mumps->myid) {
1370:       VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);
1371:       ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);
1372:     } else {
1373:       VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);
1374:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);
1375:     }
1376:     MatCreateVecs(A,NULL,&b);
1377:     VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);
1378:     ISDestroy(&is_iden);
1379:     VecDestroy(&b);
1380:     break;
1381:   }
1382:   PetscMUMPS_c(&mumps->id);
1383:   MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);

1385:   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1386:   F->ops->solve           = MatSolve_MUMPS;
1387:   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1388:   return(0);
1389: }

1391: /* Note the Petsc r permutation and factor info are ignored */
1392: PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
1393: {
1394:   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1396:   Vec            b;
1397:   IS             is_iden;
1398:   const PetscInt M = A->rmap->N;

1401:   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;

1403:   /* Set MUMPS options from the options database */
1404:   PetscSetMUMPSFromOptions(F,A);

1406:   (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);

1408:   /* analysis phase */
1409:   /*----------------*/
1410:   mumps->id.job = JOB_FACTSYMBOLIC;
1411:   mumps->id.n   = M;
1412:   switch (mumps->id.ICNTL(18)) {
1413:   case 0:  /* centralized assembled matrix input */
1414:     if (!mumps->myid) {
1415:       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1416:       if (mumps->id.ICNTL(6)>1) {
1417:         mumps->id.a = (MumpsScalar*)mumps->val;
1418:       }
1419:     }
1420:     break;
1421:   case 3:  /* distributed assembled matrix input (size>1) */
1422:     mumps->id.nz_loc = mumps->nz;
1423:     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1424:     if (mumps->id.ICNTL(6)>1) {
1425:       mumps->id.a_loc = (MumpsScalar*)mumps->val;
1426:     }
1427:     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1428:     if (!mumps->myid) {
1429:       VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);
1430:       ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);
1431:     } else {
1432:       VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);
1433:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);
1434:     }
1435:     MatCreateVecs(A,NULL,&b);
1436:     VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);
1437:     ISDestroy(&is_iden);
1438:     VecDestroy(&b);
1439:     break;
1440:   }
1441:   PetscMUMPS_c(&mumps->id);
1442:   MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);

1444:   F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
1445:   F->ops->solve                 = MatSolve_MUMPS;
1446:   F->ops->solvetranspose        = MatSolve_MUMPS;
1447:   F->ops->matsolve              = MatMatSolve_MUMPS;
1448: #if defined(PETSC_USE_COMPLEX)
1449:   F->ops->getinertia = NULL;
1450: #else
1451:   F->ops->getinertia = MatGetInertia_SBAIJMUMPS;
1452: #endif
1453:   return(0);
1454: }

1456: PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
1457: {
1458:   PetscErrorCode    ierr;
1459:   PetscBool         iascii;
1460:   PetscViewerFormat format;
1461:   Mat_MUMPS         *mumps=(Mat_MUMPS*)A->data;

1464:   /* check if matrix is mumps type */
1465:   if (A->ops->solve != MatSolve_MUMPS) return(0);

1467:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1468:   if (iascii) {
1469:     PetscViewerGetFormat(viewer,&format);
1470:     if (format == PETSC_VIEWER_ASCII_INFO) {
1471:       PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");
1472:       PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                   %d \n",mumps->id.sym);
1473:       PetscViewerASCIIPrintf(viewer,"  PAR (host participation):            %d \n",mumps->id.par);
1474:       PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):         %d \n",mumps->id.ICNTL(1));
1475:       PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg): %d \n",mumps->id.ICNTL(2));
1476:       PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):   %d \n",mumps->id.ICNTL(3));
1477:       PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):        %d \n",mumps->id.ICNTL(4));
1478:       PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):         %d \n",mumps->id.ICNTL(5));
1479:       PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):        %d \n",mumps->id.ICNTL(6));
1480:       PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (sequential matrix ordering):%d \n",mumps->id.ICNTL(7));
1481:       PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scaling strategy):        %d \n",mumps->id.ICNTL(8));
1482:       PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements):  %d \n",mumps->id.ICNTL(10));
1483:       PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):          %d \n",mumps->id.ICNTL(11));
1484:       if (mumps->id.ICNTL(11)>0) {
1485:         PetscViewerASCIIPrintf(viewer,"    RINFOG(4) (inf norm of input mat):        %g\n",mumps->id.RINFOG(4));
1486:         PetscViewerASCIIPrintf(viewer,"    RINFOG(5) (inf norm of solution):         %g\n",mumps->id.RINFOG(5));
1487:         PetscViewerASCIIPrintf(viewer,"    RINFOG(6) (inf norm of residual):         %g\n",mumps->id.RINFOG(6));
1488:         PetscViewerASCIIPrintf(viewer,"    RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8));
1489:         PetscViewerASCIIPrintf(viewer,"    RINFOG(9) (error estimate):               %g \n",mumps->id.RINFOG(9));
1490:         PetscViewerASCIIPrintf(viewer,"    RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11));
1491:       }
1492:       PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",mumps->id.ICNTL(12));
1493:       PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",mumps->id.ICNTL(13));
1494:       PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",mumps->id.ICNTL(14));
1495:       /* ICNTL(15-17) not used */
1496:       PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",mumps->id.ICNTL(18));
1497:       PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Schur complement info):                       %d \n",mumps->id.ICNTL(19));
1498:       PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",mumps->id.ICNTL(20));
1499:       PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d \n",mumps->id.ICNTL(21));
1500:       PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",mumps->id.ICNTL(22));
1501:       PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",mumps->id.ICNTL(23));

1503:       PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",mumps->id.ICNTL(24));
1504:       PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",mumps->id.ICNTL(25));
1505:       PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",mumps->id.ICNTL(26));
1506:       PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",mumps->id.ICNTL(27));
1507:       PetscViewerASCIIPrintf(viewer,"  ICNTL(28) (use parallel or sequential ordering):        %d \n",mumps->id.ICNTL(28));
1508:       PetscViewerASCIIPrintf(viewer,"  ICNTL(29) (parallel ordering):                          %d \n",mumps->id.ICNTL(29));

1510:       PetscViewerASCIIPrintf(viewer,"  ICNTL(30) (user-specified set of entries in inv(A)):    %d \n",mumps->id.ICNTL(30));
1511:       PetscViewerASCIIPrintf(viewer,"  ICNTL(31) (factors is discarded in the solve phase):    %d \n",mumps->id.ICNTL(31));
1512:       PetscViewerASCIIPrintf(viewer,"  ICNTL(33) (compute determinant):                        %d \n",mumps->id.ICNTL(33));
1513:       PetscViewerASCIIPrintf(viewer,"  ICNTL(35) (activate BLR based factorization):           %d \n",mumps->id.ICNTL(35));

1515:       PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",mumps->id.CNTL(1));
1516:       PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2));
1517:       PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",mumps->id.CNTL(3));
1518:       PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",mumps->id.CNTL(4));
1519:       PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",mumps->id.CNTL(5));
1520:       PetscViewerASCIIPrintf(viewer,"  CNTL(7) (dropping parameter for BLR):       %g \n",mumps->id.CNTL(7));

1522:       /* infomation local to each processor */
1523:       PetscViewerASCIIPrintf(viewer, "  RINFO(1) (local estimated flops for the elimination after analysis): \n");
1524:       PetscViewerASCIIPushSynchronized(viewer);
1525:       PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %g \n",mumps->myid,mumps->id.RINFO(1));
1526:       PetscViewerFlush(viewer);
1527:       PetscViewerASCIIPrintf(viewer, "  RINFO(2) (local estimated flops for the assembly after factorization): \n");
1528:       PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(2));
1529:       PetscViewerFlush(viewer);
1530:       PetscViewerASCIIPrintf(viewer, "  RINFO(3) (local estimated flops for the elimination after factorization): \n");
1531:       PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(3));
1532:       PetscViewerFlush(viewer);

1534:       PetscViewerASCIIPrintf(viewer, "  INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");
1535:       PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] %d \n",mumps->myid,mumps->id.INFO(15));
1536:       PetscViewerFlush(viewer);

1538:       PetscViewerASCIIPrintf(viewer, "  INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");
1539:       PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(16));
1540:       PetscViewerFlush(viewer);

1542:       PetscViewerASCIIPrintf(viewer, "  INFO(23) (num of pivots eliminated on this processor after factorization): \n");
1543:       PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(23));
1544:       PetscViewerFlush(viewer);

1546:       if (mumps->ninfo && mumps->ninfo <= 40){
1547:         PetscInt i;
1548:         for (i=0; i<mumps->ninfo; i++){
1549:           PetscViewerASCIIPrintf(viewer, "  INFO(%d): \n",mumps->info[i]);
1550:           PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(mumps->info[i]));
1551:           PetscViewerFlush(viewer);
1552:         }
1553:       }


1556:       PetscViewerASCIIPopSynchronized(viewer);

1558:       if (!mumps->myid) { /* information from the host */
1559:         PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",mumps->id.RINFOG(1));
1560:         PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",mumps->id.RINFOG(2));
1561:         PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",mumps->id.RINFOG(3));
1562:         PetscViewerASCIIPrintf(viewer,"  (RINFOG(12) RINFOG(13))*2^INFOG(34) (determinant): (%g,%g)*(2^%d)\n",mumps->id.RINFOG(12),mumps->id.RINFOG(13),mumps->id.INFOG(34));

1564:         PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(3));
1565:         PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(4));
1566:         PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",mumps->id.INFOG(5));
1567:         PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",mumps->id.INFOG(6));
1568:         PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively use after analysis): %d \n",mumps->id.INFOG(7));
1569:         PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",mumps->id.INFOG(8));
1570:         PetscViewerASCIIPrintf(viewer,"  INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",mumps->id.INFOG(9));
1571:         PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",mumps->id.INFOG(10));
1572:         PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",mumps->id.INFOG(11));
1573:         PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",mumps->id.INFOG(12));
1574:         PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",mumps->id.INFOG(13));
1575:         PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",mumps->id.INFOG(14));
1576:         PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",mumps->id.INFOG(15));
1577:         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",mumps->id.INFOG(16));
1578:         PetscViewerASCIIPrintf(viewer,"  INFOG(17) (estimated size of all MUMPS internal data for factorization after analysis: sum over all processors): %d \n",mumps->id.INFOG(17));
1579:         PetscViewerASCIIPrintf(viewer,"  INFOG(18) (size of all MUMPS internal data allocated during factorization: value on the most memory consuming processor): %d \n",mumps->id.INFOG(18));
1580:         PetscViewerASCIIPrintf(viewer,"  INFOG(19) (size of all MUMPS internal data allocated during factorization: sum over all processors): %d \n",mumps->id.INFOG(19));
1581:         PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",mumps->id.INFOG(20));
1582:         PetscViewerASCIIPrintf(viewer,"  INFOG(21) (size in MB of memory effectively used during factorization - value on the most memory consuming processor): %d \n",mumps->id.INFOG(21));
1583:         PetscViewerASCIIPrintf(viewer,"  INFOG(22) (size in MB of memory effectively used during factorization - sum over all processors): %d \n",mumps->id.INFOG(22));
1584:         PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",mumps->id.INFOG(23));
1585:         PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",mumps->id.INFOG(24));
1586:         PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",mumps->id.INFOG(25));
1587:         PetscViewerASCIIPrintf(viewer,"  INFOG(28) (after factorization: number of null pivots encountered): %d\n",mumps->id.INFOG(28));
1588:         PetscViewerASCIIPrintf(viewer,"  INFOG(29) (after factorization: effective number of entries in the factors (sum over all processors)): %d\n",mumps->id.INFOG(29));
1589:         PetscViewerASCIIPrintf(viewer,"  INFOG(30, 31) (after solution: size in Mbytes of memory used during solution phase): %d, %d\n",mumps->id.INFOG(30),mumps->id.INFOG(31));
1590:         PetscViewerASCIIPrintf(viewer,"  INFOG(32) (after analysis: type of analysis done): %d\n",mumps->id.INFOG(32));
1591:         PetscViewerASCIIPrintf(viewer,"  INFOG(33) (value used for ICNTL(8)): %d\n",mumps->id.INFOG(33));
1592:         PetscViewerASCIIPrintf(viewer,"  INFOG(34) (exponent of the determinant if determinant is requested): %d\n",mumps->id.INFOG(34));
1593:       }
1594:     }
1595:   }
1596:   return(0);
1597: }

1599: PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
1600: {
1601:   Mat_MUMPS *mumps =(Mat_MUMPS*)A->data;

1604:   info->block_size        = 1.0;
1605:   info->nz_allocated      = mumps->id.INFOG(20);
1606:   info->nz_used           = mumps->id.INFOG(20);
1607:   info->nz_unneeded       = 0.0;
1608:   info->assemblies        = 0.0;
1609:   info->mallocs           = 0.0;
1610:   info->memory            = 0.0;
1611:   info->fill_ratio_given  = 0;
1612:   info->fill_ratio_needed = 0;
1613:   info->factor_mallocs    = 0;
1614:   return(0);
1615: }

1617: /* -------------------------------------------------------------------------------------------*/
1618: PetscErrorCode MatFactorSetSchurIS_MUMPS(Mat F, IS is)
1619: {
1620:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1621:   const PetscInt *idxs;
1622:   PetscInt       size,i;

1626:   ISGetLocalSize(is,&size);
1627:   if (mumps->size > 1) {
1628:     PetscBool ls,gs;

1630:     ls   = mumps->myid ? (size ? PETSC_FALSE : PETSC_TRUE) : PETSC_TRUE;
1631:     MPI_Allreduce(&ls,&gs,1,MPIU_BOOL,MPI_LAND,mumps->comm_mumps);
1632:     if (!gs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MUMPS distributed parallel Schur complements not yet supported from PETSc\n");
1633:   }
1634:   if (mumps->id.size_schur != size) {
1635:     PetscFree2(mumps->id.listvar_schur,mumps->id.schur);
1636:     mumps->id.size_schur = size;
1637:     mumps->id.schur_lld  = size;
1638:     PetscMalloc2(size,&mumps->id.listvar_schur,size*size,&mumps->id.schur);
1639:   }

1641:   /* Schur complement matrix */
1642:   MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.size_schur,(PetscScalar*)mumps->id.schur,&F->schur);
1643:   if (mumps->sym == 1) {
1644:     MatSetOption(F->schur,MAT_SPD,PETSC_TRUE);
1645:   }

1647:   /* MUMPS expects Fortran style indices */
1648:   ISGetIndices(is,&idxs);
1649:   PetscMemcpy(mumps->id.listvar_schur,idxs,size*sizeof(PetscInt));
1650:   for (i=0;i<size;i++) mumps->id.listvar_schur[i]++;
1651:   ISRestoreIndices(is,&idxs);
1652:   if (mumps->size > 1) {
1653:     mumps->id.ICNTL(19) = 1; /* MUMPS returns Schur centralized on the host */
1654:   } else {
1655:     if (F->factortype == MAT_FACTOR_LU) {
1656:       mumps->id.ICNTL(19) = 3; /* MUMPS returns full matrix */
1657:     } else {
1658:       mumps->id.ICNTL(19) = 2; /* MUMPS returns lower triangular part */
1659:     }
1660:   }
1661:   /* set a special value of ICNTL (not handled my MUMPS) to be used in the solve phase by PETSc */
1662:   mumps->id.ICNTL(26) = -1;
1663:   return(0);
1664: }

1666: /* -------------------------------------------------------------------------------------------*/
1667: PetscErrorCode MatFactorCreateSchurComplement_MUMPS(Mat F,Mat* S)
1668: {
1669:   Mat            St;
1670:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1671:   PetscScalar    *array;
1672: #if defined(PETSC_USE_COMPLEX)
1673:   PetscScalar    im = PetscSqrtScalar((PetscScalar)-1.0);
1674: #endif

1678:   if (!mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it");
1679:   MatCreate(PETSC_COMM_SELF,&St);
1680:   MatSetSizes(St,PETSC_DECIDE,PETSC_DECIDE,mumps->id.size_schur,mumps->id.size_schur);
1681:   MatSetType(St,MATDENSE);
1682:   MatSetUp(St);
1683:   MatDenseGetArray(St,&array);
1684:   if (!mumps->sym) { /* MUMPS always return a full matrix */
1685:     if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
1686:       PetscInt i,j,N=mumps->id.size_schur;
1687:       for (i=0;i<N;i++) {
1688:         for (j=0;j<N;j++) {
1689: #if !defined(PETSC_USE_COMPLEX)
1690:           PetscScalar val = mumps->id.schur[i*N+j];
1691: #else
1692:           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
1693: #endif
1694:           array[j*N+i] = val;
1695:         }
1696:       }
1697:     } else { /* stored by columns */
1698:       PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));
1699:     }
1700:   } else { /* either full or lower-triangular (not packed) */
1701:     if (mumps->id.ICNTL(19) == 2) { /* lower triangular stored by columns */
1702:       PetscInt i,j,N=mumps->id.size_schur;
1703:       for (i=0;i<N;i++) {
1704:         for (j=i;j<N;j++) {
1705: #if !defined(PETSC_USE_COMPLEX)
1706:           PetscScalar val = mumps->id.schur[i*N+j];
1707: #else
1708:           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
1709: #endif
1710:           array[i*N+j] = val;
1711:           array[j*N+i] = val;
1712:         }
1713:       }
1714:     } else if (mumps->id.ICNTL(19) == 3) { /* full matrix */
1715:       PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));
1716:     } else { /* ICNTL(19) == 1 lower triangular stored by rows */
1717:       PetscInt i,j,N=mumps->id.size_schur;
1718:       for (i=0;i<N;i++) {
1719:         for (j=0;j<i+1;j++) {
1720: #if !defined(PETSC_USE_COMPLEX)
1721:           PetscScalar val = mumps->id.schur[i*N+j];
1722: #else
1723:           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
1724: #endif
1725:           array[i*N+j] = val;
1726:           array[j*N+i] = val;
1727:         }
1728:       }
1729:     }
1730:   }
1731:   MatDenseRestoreArray(St,&array);
1732:   *S   = St;
1733:   return(0);
1734: }

1736: /* -------------------------------------------------------------------------------------------*/
1737: PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival)
1738: {
1739:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

1742:   mumps->id.ICNTL(icntl) = ival;
1743:   return(0);
1744: }

1746: PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt *ival)
1747: {
1748:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

1751:   *ival = mumps->id.ICNTL(icntl);
1752:   return(0);
1753: }

1755: /*@
1756:   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()

1758:    Logically Collective on Mat

1760:    Input Parameters:
1761: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1762: .  icntl - index of MUMPS parameter array ICNTL()
1763: -  ival - value of MUMPS ICNTL(icntl)

1765:   Options Database:
1766: .   -mat_mumps_icntl_<icntl> <ival>

1768:    Level: beginner

1770:    References:
1771: .     MUMPS Users' Guide

1773: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
1774:  @*/
1775: PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival)
1776: {
1778: 
1781:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
1784:   PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));
1785:   return(0);
1786: }

1788: /*@
1789:   MatMumpsGetIcntl - Get MUMPS parameter ICNTL()

1791:    Logically Collective on Mat

1793:    Input Parameters:
1794: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1795: -  icntl - index of MUMPS parameter array ICNTL()

1797:   Output Parameter:
1798: .  ival - value of MUMPS ICNTL(icntl)

1800:    Level: beginner

1802:    References:
1803: .     MUMPS Users' Guide

1805: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
1806: @*/
1807: PetscErrorCode MatMumpsGetIcntl(Mat F,PetscInt icntl,PetscInt *ival)
1808: {

1813:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
1816:   PetscUseMethod(F,"MatMumpsGetIcntl_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));
1817:   return(0);
1818: }

1820: /* -------------------------------------------------------------------------------------------*/
1821: PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal val)
1822: {
1823:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

1826:   mumps->id.CNTL(icntl) = val;
1827:   return(0);
1828: }

1830: PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal *val)
1831: {
1832:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

1835:   *val = mumps->id.CNTL(icntl);
1836:   return(0);
1837: }

1839: /*@
1840:   MatMumpsSetCntl - Set MUMPS parameter CNTL()

1842:    Logically Collective on Mat

1844:    Input Parameters:
1845: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1846: .  icntl - index of MUMPS parameter array CNTL()
1847: -  val - value of MUMPS CNTL(icntl)

1849:   Options Database:
1850: .   -mat_mumps_cntl_<icntl> <val>

1852:    Level: beginner

1854:    References:
1855: .     MUMPS Users' Guide

1857: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
1858: @*/
1859: PetscErrorCode MatMumpsSetCntl(Mat F,PetscInt icntl,PetscReal val)
1860: {

1865:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
1868:   PetscTryMethod(F,"MatMumpsSetCntl_C",(Mat,PetscInt,PetscReal),(F,icntl,val));
1869:   return(0);
1870: }

1872: /*@
1873:   MatMumpsGetCntl - Get MUMPS parameter CNTL()

1875:    Logically Collective on Mat

1877:    Input Parameters:
1878: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1879: -  icntl - index of MUMPS parameter array CNTL()

1881:   Output Parameter:
1882: .  val - value of MUMPS CNTL(icntl)

1884:    Level: beginner

1886:    References:
1887: .      MUMPS Users' Guide

1889: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
1890: @*/
1891: PetscErrorCode MatMumpsGetCntl(Mat F,PetscInt icntl,PetscReal *val)
1892: {

1897:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
1900:   PetscUseMethod(F,"MatMumpsGetCntl_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));
1901:   return(0);
1902: }

1904: PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F,PetscInt icntl,PetscInt *info)
1905: {
1906:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

1909:   *info = mumps->id.INFO(icntl);
1910:   return(0);
1911: }

1913: PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F,PetscInt icntl,PetscInt *infog)
1914: {
1915:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

1918:   *infog = mumps->id.INFOG(icntl);
1919:   return(0);
1920: }

1922: PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfo)
1923: {
1924:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

1927:   *rinfo = mumps->id.RINFO(icntl);
1928:   return(0);
1929: }

1931: PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfog)
1932: {
1933:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

1936:   *rinfog = mumps->id.RINFOG(icntl);
1937:   return(0);
1938: }

1940: PetscErrorCode MatMumpsGetInverse_MUMPS(Mat F,Mat spRHS)
1941: {
1943:   Mat            Bt = NULL;
1944:   PetscBool      flgT;
1945:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1946:   PetscBool      done;
1947:   PetscScalar    *aa;
1948:   PetscInt       spnr,*ia,*ja;

1951:   if (!mumps->myid) {
1953:     PetscObjectTypeCompare((PetscObject)spRHS,MATTRANSPOSEMAT,&flgT);
1954:     if (flgT) {
1955:       MatTransposeGetMat(spRHS,&Bt);
1956:     } else {
1957:       SETERRQ(PetscObjectComm((PetscObject)spRHS),PETSC_ERR_ARG_WRONG,"Matrix spRHS must be type MATTRANSPOSEMAT matrix");
1958:     }
1959:   }

1961:   MatMumpsSetIcntl(F,30,1);

1963:   if (!mumps->myid) {
1964:     MatSeqAIJGetArray(Bt,&aa);
1965:     MatGetRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&done);
1966:     if (!done) SETERRQ(PetscObjectComm((PetscObject)Bt),PETSC_ERR_ARG_WRONG,"Cannot get IJ structure");

1968:     mumps->id.irhs_ptr    = ia;
1969:     mumps->id.irhs_sparse = ja;
1970:     mumps->id.nz_rhs      = ia[spnr] - 1;
1971:     mumps->id.rhs_sparse  = (MumpsScalar*)aa;
1972:   } else {
1973:     mumps->id.irhs_ptr    = NULL;
1974:     mumps->id.irhs_sparse = NULL;
1975:     mumps->id.nz_rhs      = 0;
1976:     mumps->id.rhs_sparse  = NULL;
1977:   }
1978:   mumps->id.ICNTL(20)   = 1; /* rhs is sparse */
1979:   mumps->id.ICNTL(21)   = 0; /* solution is in assembled centralized format */

1981:   /* solve phase */
1982:   /*-------------*/
1983:   mumps->id.job = JOB_SOLVE;
1984:   PetscMUMPS_c(&mumps->id);
1985:   if (mumps->id.INFOG(1) < 0)
1986:     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));

1988:   if (!mumps->myid) {
1989:     MatSeqAIJRestoreArray(Bt,&aa);
1990:     MatRestoreRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&done);
1991:   }
1992:   return(0);
1993: }

1995: /*@
1996:   MatMumpsGetInverse - Get user-specified set of entries in inverse of A

1998:    Logically Collective on Mat

2000:    Input Parameters:
2001: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2002: -  spRHS - sequential sparse matrix in MATTRANSPOSEMAT format holding specified indices in processor[0]

2004:   Output Parameter:
2005: . spRHS - requested entries of inverse of A

2007:    Level: beginner

2009:    References:
2010: .      MUMPS Users' Guide

2012: .seealso: MatGetFactor(), MatCreateTranspose()
2013: @*/
2014: PetscErrorCode MatMumpsGetInverse(Mat F,Mat spRHS)
2015: {

2020:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2021:   PetscUseMethod(F,"MatMumpsGetInverse_C",(Mat,Mat),(F,spRHS));
2022:   return(0);
2023: }

2025: /*@
2026:   MatMumpsGetInfo - Get MUMPS parameter INFO()

2028:    Logically Collective on Mat

2030:    Input Parameters:
2031: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2032: -  icntl - index of MUMPS parameter array INFO()

2034:   Output Parameter:
2035: .  ival - value of MUMPS INFO(icntl)

2037:    Level: beginner

2039:    References:
2040: .      MUMPS Users' Guide

2042: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2043: @*/
2044: PetscErrorCode MatMumpsGetInfo(Mat F,PetscInt icntl,PetscInt *ival)
2045: {

2050:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2052:   PetscUseMethod(F,"MatMumpsGetInfo_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));
2053:   return(0);
2054: }

2056: /*@
2057:   MatMumpsGetInfog - Get MUMPS parameter INFOG()

2059:    Logically Collective on Mat

2061:    Input Parameters:
2062: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2063: -  icntl - index of MUMPS parameter array INFOG()

2065:   Output Parameter:
2066: .  ival - value of MUMPS INFOG(icntl)

2068:    Level: beginner

2070:    References:
2071: .      MUMPS Users' Guide

2073: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2074: @*/
2075: PetscErrorCode MatMumpsGetInfog(Mat F,PetscInt icntl,PetscInt *ival)
2076: {

2081:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2083:   PetscUseMethod(F,"MatMumpsGetInfog_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));
2084:   return(0);
2085: }

2087: /*@
2088:   MatMumpsGetRinfo - Get MUMPS parameter RINFO()

2090:    Logically Collective on Mat

2092:    Input Parameters:
2093: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2094: -  icntl - index of MUMPS parameter array RINFO()

2096:   Output Parameter:
2097: .  val - value of MUMPS RINFO(icntl)

2099:    Level: beginner

2101:    References:
2102: .       MUMPS Users' Guide

2104: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2105: @*/
2106: PetscErrorCode MatMumpsGetRinfo(Mat F,PetscInt icntl,PetscReal *val)
2107: {

2112:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2114:   PetscUseMethod(F,"MatMumpsGetRinfo_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));
2115:   return(0);
2116: }

2118: /*@
2119:   MatMumpsGetRinfog - Get MUMPS parameter RINFOG()

2121:    Logically Collective on Mat

2123:    Input Parameters:
2124: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2125: -  icntl - index of MUMPS parameter array RINFOG()

2127:   Output Parameter:
2128: .  val - value of MUMPS RINFOG(icntl)

2130:    Level: beginner

2132:    References:
2133: .      MUMPS Users' Guide

2135: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2136: @*/
2137: PetscErrorCode MatMumpsGetRinfog(Mat F,PetscInt icntl,PetscReal *val)
2138: {

2143:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2145:   PetscUseMethod(F,"MatMumpsGetRinfog_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));
2146:   return(0);
2147: }

2149: /*MC
2150:   MATSOLVERMUMPS -  A matrix type providing direct solvers (LU and Cholesky) for
2151:   distributed and sequential matrices via the external package MUMPS.

2153:   Works with MATAIJ and MATSBAIJ matrices

2155:   Use ./configure --download-mumps --download-scalapack --download-parmetis --download-metis --download-ptscotch  to have PETSc installed with MUMPS

2157:   Use -pc_type cholesky or lu -pc_factor_mat_solver_type mumps to use this direct solver

2159:   Options Database Keys:
2160: +  -mat_mumps_icntl_1 - ICNTL(1): output stream for error messages 
2161: .  -mat_mumps_icntl_2 - ICNTL(2): output stream for diagnostic printing, statistics, and warning 
2162: .  -mat_mumps_icntl_3 -  ICNTL(3): output stream for global information, collected on the host 
2163: .  -mat_mumps_icntl_4 -  ICNTL(4): level of printing (0 to 4) 
2164: .  -mat_mumps_icntl_6 - ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7) 
2165: .  -mat_mumps_icntl_7 - ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis 
2166: .  -mat_mumps_icntl_8  - ICNTL(8): scaling strategy (-2 to 8 or 77) 
2167: .  -mat_mumps_icntl_10  - ICNTL(10): max num of refinements 
2168: .  -mat_mumps_icntl_11  - ICNTL(11): statistics related to an error analysis (via -ksp_view) 
2169: .  -mat_mumps_icntl_12  - ICNTL(12): an ordering strategy for symmetric matrices (0 to 3) 
2170: .  -mat_mumps_icntl_13  - ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting 
2171: .  -mat_mumps_icntl_14  - ICNTL(14): percentage increase in the estimated working space 
2172: .  -mat_mumps_icntl_19  - ICNTL(19): computes the Schur complement 
2173: .  -mat_mumps_icntl_22  - ICNTL(22): in-core/out-of-core factorization and solve (0 or 1) 
2174: .  -mat_mumps_icntl_23  - ICNTL(23): max size of the working memory (MB) that can allocate per processor 
2175: .  -mat_mumps_icntl_24  - ICNTL(24): detection of null pivot rows (0 or 1) 
2176: .  -mat_mumps_icntl_25  - ICNTL(25): compute a solution of a deficient matrix and a null space basis 
2177: .  -mat_mumps_icntl_26  - ICNTL(26): drives the solution phase if a Schur complement matrix 
2178: .  -mat_mumps_icntl_28  - ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering 
2179: .  -mat_mumps_icntl_29 - ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis 
2180: .  -mat_mumps_icntl_30 - ICNTL(30): compute user-specified set of entries in inv(A) 
2181: .  -mat_mumps_icntl_31 - ICNTL(31): indicates which factors may be discarded during factorization 
2182: .  -mat_mumps_icntl_33 - ICNTL(33): compute determinant 
2183: .  -mat_mumps_cntl_1  - CNTL(1): relative pivoting threshold 
2184: .  -mat_mumps_cntl_2  -  CNTL(2): stopping criterion of refinement 
2185: .  -mat_mumps_cntl_3 - CNTL(3): absolute pivoting threshold 
2186: .  -mat_mumps_cntl_4 - CNTL(4): value for static pivoting 
2187: -  -mat_mumps_cntl_5 - CNTL(5): fixation for null pivots 

2189:   Level: beginner

2191:     Notes: When a MUMPS factorization fails inside a KSP solve, for example with a KSP_DIVERGED_PCSETUP_FAILED, one can find the MUMPS information about the failure by calling 
2192: $          KSPGetPC(ksp,&pc);
2193: $          PCFactorGetMatrix(pc,&mat);
2194: $          MatMumpsGetInfo(mat,....);
2195: $          MatMumpsGetInfog(mat,....); etc.
2196:            Or you can run with -ksp_error_if_not_converged and the program will be stopped and the information printed in the error message.

2198: .seealso: PCFactorSetMatSolverType(), MatSolverType, MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog(), KSPGetPC(), PCGetFactor(), PCFactorGetMatrix()

2200: M*/

2202: static PetscErrorCode MatFactorGetSolverType_mumps(Mat A,MatSolverType *type)
2203: {
2205:   *type = MATSOLVERMUMPS;
2206:   return(0);
2207: }

2209: /* MatGetFactor for Seq and MPI AIJ matrices */
2210: static PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
2211: {
2212:   Mat            B;
2214:   Mat_MUMPS      *mumps;
2215:   PetscBool      isSeqAIJ;

2218:   /* Create the factorization matrix */
2219:   PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);
2220:   MatCreate(PetscObjectComm((PetscObject)A),&B);
2221:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
2222:   PetscStrallocpy("mumps",&((PetscObject)B)->type_name);
2223:   MatSetUp(B);

2225:   PetscNewLog(B,&mumps);

2227:   B->ops->view        = MatView_MUMPS;
2228:   B->ops->getinfo     = MatGetInfo_MUMPS;

2230:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);
2231:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);
2232:   PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);
2233:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);
2234:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);
2235:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);
2236:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);
2237:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);
2238:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);
2239:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);
2240:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);
2241:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverse_C",MatMumpsGetInverse_MUMPS);

2243:   if (ftype == MAT_FACTOR_LU) {
2244:     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
2245:     B->factortype            = MAT_FACTOR_LU;
2246:     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
2247:     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
2248:     mumps->sym = 0;
2249:   } else {
2250:     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2251:     B->factortype                  = MAT_FACTOR_CHOLESKY;
2252:     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
2253:     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
2254: #if defined(PETSC_USE_COMPLEX)
2255:     mumps->sym = 2;
2256: #else
2257:     if (A->spd_set && A->spd) mumps->sym = 1;
2258:     else                      mumps->sym = 2;
2259: #endif
2260:   }

2262:   /* set solvertype */
2263:   PetscFree(B->solvertype);
2264:   PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);

2266:   mumps->isAIJ    = PETSC_TRUE;
2267:   B->ops->destroy = MatDestroy_MUMPS;
2268:   B->data        = (void*)mumps;

2270:   PetscInitializeMUMPS(A,mumps);

2272:   *F = B;
2273:   return(0);
2274: }

2276: /* MatGetFactor for Seq and MPI SBAIJ matrices */
2277: static PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
2278: {
2279:   Mat            B;
2281:   Mat_MUMPS      *mumps;
2282:   PetscBool      isSeqSBAIJ;

2285:   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
2286:   if (A->rmap->bs > 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with block size > 1 with MUMPS Cholesky, use AIJ matrix instead");
2287:   PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);
2288:   /* Create the factorization matrix */
2289:   MatCreate(PetscObjectComm((PetscObject)A),&B);
2290:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
2291:   PetscStrallocpy("mumps",&((PetscObject)B)->type_name);
2292:   MatSetUp(B);

2294:   PetscNewLog(B,&mumps);
2295:   if (isSeqSBAIJ) {
2296:     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
2297:   } else {
2298:     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
2299:   }

2301:   B->ops->getinfo                = MatGetInfo_External;
2302:   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2303:   B->ops->view                   = MatView_MUMPS;

2305:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);
2306:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);
2307:   PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);
2308:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);
2309:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);
2310:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);
2311:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);
2312:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);
2313:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);
2314:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);
2315:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);
2316:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverse_C",MatMumpsGetInverse_MUMPS);

2318:   B->factortype = MAT_FACTOR_CHOLESKY;
2319: #if defined(PETSC_USE_COMPLEX)
2320:   mumps->sym = 2;
2321: #else
2322:   if (A->spd_set && A->spd) mumps->sym = 1;
2323:   else                      mumps->sym = 2;
2324: #endif

2326:   /* set solvertype */
2327:   PetscFree(B->solvertype);
2328:   PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);

2330:   mumps->isAIJ    = PETSC_FALSE;
2331:   B->ops->destroy = MatDestroy_MUMPS;
2332:   B->data        = (void*)mumps;

2334:   PetscInitializeMUMPS(A,mumps);

2336:   *F = B;
2337:   return(0);
2338: }

2340: static PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
2341: {
2342:   Mat            B;
2344:   Mat_MUMPS      *mumps;
2345:   PetscBool      isSeqBAIJ;

2348:   /* Create the factorization matrix */
2349:   PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);
2350:   MatCreate(PetscObjectComm((PetscObject)A),&B);
2351:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
2352:   PetscStrallocpy("mumps",&((PetscObject)B)->type_name);
2353:   MatSetUp(B);

2355:   PetscNewLog(B,&mumps);
2356:   if (ftype == MAT_FACTOR_LU) {
2357:     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
2358:     B->factortype            = MAT_FACTOR_LU;
2359:     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
2360:     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
2361:     mumps->sym = 0;
2362:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n");

2364:   B->ops->getinfo     = MatGetInfo_External;
2365:   B->ops->view        = MatView_MUMPS;

2367:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);
2368:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);
2369:   PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);
2370:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);
2371:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);
2372:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);
2373:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);
2374:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);
2375:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);
2376:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);
2377:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);
2378:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverse_C",MatMumpsGetInverse_MUMPS);

2380:   /* set solvertype */
2381:   PetscFree(B->solvertype);
2382:   PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);

2384:   mumps->isAIJ    = PETSC_TRUE;
2385:   B->ops->destroy = MatDestroy_MUMPS;
2386:   B->data        = (void*)mumps;

2388:   PetscInitializeMUMPS(A,mumps);

2390:   *F = B;
2391:   return(0);
2392: }

2394: PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_MUMPS(void)
2395: {

2399:   MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);
2400:   MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);
2401:   MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);
2402:   MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);
2403:   MatSolverTypeRegister(MATSOLVERMUMPS,MATMPISBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);
2404:   MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);
2405:   MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);
2406:   MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);
2407:   MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);
2408:   MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);
2409:   return(0);
2410: }