Actual source code: mumps.c

petsc-3.9.4 2018-09-11
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>
  8:  #include <../src/mat/impls/sell/mpi/mpisell.h>

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

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

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

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

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

 80:   MatStructure matstruc;
 81:   PetscMPIInt  myid,size;
 82:   PetscInt     *irn,*jcn,nz,sym;
 83:   PetscScalar  *val;
 84:   MPI_Comm     comm_mumps;
 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_seqsell_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
246: {
247:   Mat_SeqSELL *a=(Mat_SeqSELL*)A->data;
248:   PetscInt    *ptr;

251:   *v = a->val;
252:   if (reuse == MAT_INITIAL_MATRIX) {
253:     PetscInt       nz,i,j,row;

256:     nz   = a->sliidx[a->totalslices];
257:     *nnz = nz;
258:     PetscMalloc1(2*nz, &ptr);
259:     *r   = ptr;
260:     *c   = ptr + nz;

262:     for (i=0; i<a->totalslices; i++) {
263:       for (j=a->sliidx[i],row=0; j<a->sliidx[i+1]; j++,row=((row+1)&0x07)) {
264:         *ptr++ = 8*i + row + shift;
265:       }
266:     }
267:     for (i=0;i<nz;i++) *ptr++ = a->colidx[i] + shift;
268:   }
269:   return(0);
270: }

272: PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
273: {
274:   Mat_SeqBAIJ    *aa=(Mat_SeqBAIJ*)A->data;
275:   const PetscInt *ai,*aj,*ajj,bs2 = aa->bs2;
276:   PetscInt       bs,M,nz,idx=0,rnz,i,j,k,m;
278:   PetscInt       *row,*col;

281:   MatGetBlockSize(A,&bs);
282:   M = A->rmap->N/bs;
283:   *v = aa->a;
284:   if (reuse == MAT_INITIAL_MATRIX) {
285:     ai   = aa->i; aj = aa->j;
286:     nz   = bs2*aa->nz;
287:     *nnz = nz;
288:     PetscMalloc1(2*nz, &row);
289:     col  = row + nz;

291:     for (i=0; i<M; i++) {
292:       ajj = aj + ai[i];
293:       rnz = ai[i+1] - ai[i];
294:       for (k=0; k<rnz; k++) {
295:         for (j=0; j<bs; j++) {
296:           for (m=0; m<bs; m++) {
297:             row[idx]   = i*bs + m + shift;
298:             col[idx++] = bs*(ajj[k]) + j + shift;
299:           }
300:         }
301:       }
302:     }
303:     *r = row; *c = col;
304:   }
305:   return(0);
306: }

308: PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
309: {
310:   const PetscInt *ai, *aj,*ajj,M=A->rmap->n;
311:   PetscInt       nz,rnz,i,j;
313:   PetscInt       *row,*col;
314:   Mat_SeqSBAIJ   *aa=(Mat_SeqSBAIJ*)A->data;

317:   *v = aa->a;
318:   if (reuse == MAT_INITIAL_MATRIX) {
319:     nz   = aa->nz;
320:     ai   = aa->i;
321:     aj   = aa->j;
322:     *v   = aa->a;
323:     *nnz = nz;
324:     PetscMalloc1(2*nz, &row);
325:     col  = row + nz;

327:     nz = 0;
328:     for (i=0; i<M; i++) {
329:       rnz = ai[i+1] - ai[i];
330:       ajj = aj + ai[i];
331:       for (j=0; j<rnz; j++) {
332:         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
333:       }
334:     }
335:     *r = row; *c = col;
336:   }
337:   return(0);
338: }

340: PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
341: {
342:   const PetscInt    *ai,*aj,*ajj,*adiag,M=A->rmap->n;
343:   PetscInt          nz,rnz,i,j;
344:   const PetscScalar *av,*v1;
345:   PetscScalar       *val;
346:   PetscErrorCode    ierr;
347:   PetscInt          *row,*col;
348:   Mat_SeqAIJ        *aa=(Mat_SeqAIJ*)A->data;
349:   PetscBool         missing;

352:   ai    = aa->i; aj = aa->j; av = aa->a;
353:   adiag = aa->diag;
354:   MatMissingDiagonal_SeqAIJ(A,&missing,&i);
355:   if (reuse == MAT_INITIAL_MATRIX) {
356:     /* count nz in the upper triangular part of A */
357:     nz = 0;
358:     if (missing) {
359:       for (i=0; i<M; i++) {
360:         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
361:           for (j=ai[i];j<ai[i+1];j++) {
362:             if (aj[j] < i) continue;
363:             nz++;
364:           }
365:         } else {
366:           nz += ai[i+1] - adiag[i];
367:         }
368:       }
369:     } else {
370:       for (i=0; i<M; i++) nz += ai[i+1] - adiag[i];
371:     }
372:     *nnz = nz;

374:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
375:     col  = row + nz;
376:     val  = (PetscScalar*)(col + nz);

378:     nz = 0;
379:     if (missing) {
380:       for (i=0; i<M; i++) {
381:         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
382:           for (j=ai[i];j<ai[i+1];j++) {
383:             if (aj[j] < i) continue;
384:             row[nz] = i+shift;
385:             col[nz] = aj[j]+shift;
386:             val[nz] = av[j];
387:             nz++;
388:           }
389:         } else {
390:           rnz = ai[i+1] - adiag[i];
391:           ajj = aj + adiag[i];
392:           v1  = av + adiag[i];
393:           for (j=0; j<rnz; j++) {
394:             row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
395:           }
396:         }
397:       }
398:     } else {
399:       for (i=0; i<M; i++) {
400:         rnz = ai[i+1] - adiag[i];
401:         ajj = aj + adiag[i];
402:         v1  = av + adiag[i];
403:         for (j=0; j<rnz; j++) {
404:           row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
405:         }
406:       }
407:     }
408:     *r = row; *c = col; *v = val;
409:   } else {
410:     nz = 0; val = *v;
411:     if (missing) {
412:       for (i=0; i <M; i++) {
413:         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
414:           for (j=ai[i];j<ai[i+1];j++) {
415:             if (aj[j] < i) continue;
416:             val[nz++] = av[j];
417:           }
418:         } else {
419:           rnz = ai[i+1] - adiag[i];
420:           v1  = av + adiag[i];
421:           for (j=0; j<rnz; j++) {
422:             val[nz++] = v1[j];
423:           }
424:         }
425:       }
426:     } else {
427:       for (i=0; i <M; i++) {
428:         rnz = ai[i+1] - adiag[i];
429:         v1  = av + adiag[i];
430:         for (j=0; j<rnz; j++) {
431:           val[nz++] = v1[j];
432:         }
433:       }
434:     }
435:   }
436:   return(0);
437: }

439: PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
440: {
441:   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
442:   PetscErrorCode    ierr;
443:   PetscInt          rstart,nz,i,j,jj,irow,countA,countB;
444:   PetscInt          *row,*col;
445:   const PetscScalar *av, *bv,*v1,*v2;
446:   PetscScalar       *val;
447:   Mat_MPISBAIJ      *mat = (Mat_MPISBAIJ*)A->data;
448:   Mat_SeqSBAIJ      *aa  = (Mat_SeqSBAIJ*)(mat->A)->data;
449:   Mat_SeqBAIJ       *bb  = (Mat_SeqBAIJ*)(mat->B)->data;

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

455:   garray = mat->garray;

457:   if (reuse == MAT_INITIAL_MATRIX) {
458:     nz   = aa->nz + bb->nz;
459:     *nnz = nz;
460:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
461:     col  = row + nz;
462:     val  = (PetscScalar*)(col + nz);

464:     *r = row; *c = col; *v = val;
465:   } else {
466:     row = *r; col = *c; val = *v;
467:   }

469:   jj = 0; irow = rstart;
470:   for (i=0; i<m; i++) {
471:     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
472:     countA = ai[i+1] - ai[i];
473:     countB = bi[i+1] - bi[i];
474:     bjj    = bj + bi[i];
475:     v1     = av + ai[i];
476:     v2     = bv + bi[i];

478:     /* A-part */
479:     for (j=0; j<countA; j++) {
480:       if (reuse == MAT_INITIAL_MATRIX) {
481:         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
482:       }
483:       val[jj++] = v1[j];
484:     }

486:     /* B-part */
487:     for (j=0; j < countB; j++) {
488:       if (reuse == MAT_INITIAL_MATRIX) {
489:         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
490:       }
491:       val[jj++] = v2[j];
492:     }
493:     irow++;
494:   }
495:   return(0);
496: }

498: PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
499: {
500:   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
501:   PetscErrorCode    ierr;
502:   PetscInt          rstart,nz,i,j,jj,irow,countA,countB;
503:   PetscInt          *row,*col;
504:   const PetscScalar *av, *bv,*v1,*v2;
505:   PetscScalar       *val;
506:   Mat_MPIAIJ        *mat = (Mat_MPIAIJ*)A->data;
507:   Mat_SeqAIJ        *aa  = (Mat_SeqAIJ*)(mat->A)->data;
508:   Mat_SeqAIJ        *bb  = (Mat_SeqAIJ*)(mat->B)->data;

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

514:   garray = mat->garray;

516:   if (reuse == MAT_INITIAL_MATRIX) {
517:     nz   = aa->nz + bb->nz;
518:     *nnz = nz;
519:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
520:     col  = row + nz;
521:     val  = (PetscScalar*)(col + nz);

523:     *r = row; *c = col; *v = val;
524:   } else {
525:     row = *r; col = *c; val = *v;
526:   }

528:   jj = 0; irow = rstart;
529:   for (i=0; i<m; i++) {
530:     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
531:     countA = ai[i+1] - ai[i];
532:     countB = bi[i+1] - bi[i];
533:     bjj    = bj + bi[i];
534:     v1     = av + ai[i];
535:     v2     = bv + bi[i];

537:     /* A-part */
538:     for (j=0; j<countA; j++) {
539:       if (reuse == MAT_INITIAL_MATRIX) {
540:         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
541:       }
542:       val[jj++] = v1[j];
543:     }

545:     /* B-part */
546:     for (j=0; j < countB; j++) {
547:       if (reuse == MAT_INITIAL_MATRIX) {
548:         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
549:       }
550:       val[jj++] = v2[j];
551:     }
552:     irow++;
553:   }
554:   return(0);
555: }

557: PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
558: {
559:   Mat_MPIBAIJ       *mat    = (Mat_MPIBAIJ*)A->data;
560:   Mat_SeqBAIJ       *aa     = (Mat_SeqBAIJ*)(mat->A)->data;
561:   Mat_SeqBAIJ       *bb     = (Mat_SeqBAIJ*)(mat->B)->data;
562:   const PetscInt    *ai     = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj;
563:   const PetscInt    *garray = mat->garray,mbs=mat->mbs,rstart=A->rmap->rstart;
564:   const PetscInt    bs2=mat->bs2;
565:   PetscErrorCode    ierr;
566:   PetscInt          bs,nz,i,j,k,n,jj,irow,countA,countB,idx;
567:   PetscInt          *row,*col;
568:   const PetscScalar *av=aa->a, *bv=bb->a,*v1,*v2;
569:   PetscScalar       *val;

572:   MatGetBlockSize(A,&bs);
573:   if (reuse == MAT_INITIAL_MATRIX) {
574:     nz   = bs2*(aa->nz + bb->nz);
575:     *nnz = nz;
576:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
577:     col  = row + nz;
578:     val  = (PetscScalar*)(col + nz);

580:     *r = row; *c = col; *v = val;
581:   } else {
582:     row = *r; col = *c; val = *v;
583:   }

585:   jj = 0; irow = rstart;
586:   for (i=0; i<mbs; i++) {
587:     countA = ai[i+1] - ai[i];
588:     countB = bi[i+1] - bi[i];
589:     ajj    = aj + ai[i];
590:     bjj    = bj + bi[i];
591:     v1     = av + bs2*ai[i];
592:     v2     = bv + bs2*bi[i];

594:     idx = 0;
595:     /* A-part */
596:     for (k=0; k<countA; k++) {
597:       for (j=0; j<bs; j++) {
598:         for (n=0; n<bs; n++) {
599:           if (reuse == MAT_INITIAL_MATRIX) {
600:             row[jj] = irow + n + shift;
601:             col[jj] = rstart + bs*ajj[k] + j + shift;
602:           }
603:           val[jj++] = v1[idx++];
604:         }
605:       }
606:     }

608:     idx = 0;
609:     /* B-part */
610:     for (k=0; k<countB; k++) {
611:       for (j=0; j<bs; j++) {
612:         for (n=0; n<bs; n++) {
613:           if (reuse == MAT_INITIAL_MATRIX) {
614:             row[jj] = irow + n + shift;
615:             col[jj] = bs*garray[bjj[k]] + j + shift;
616:           }
617:           val[jj++] = v2[idx++];
618:         }
619:       }
620:     }
621:     irow += bs;
622:   }
623:   return(0);
624: }

626: PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
627: {
628:   const PetscInt    *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
629:   PetscErrorCode    ierr;
630:   PetscInt          rstart,nz,nza,nzb,i,j,jj,irow,countA,countB;
631:   PetscInt          *row,*col;
632:   const PetscScalar *av, *bv,*v1,*v2;
633:   PetscScalar       *val;
634:   Mat_MPIAIJ        *mat =  (Mat_MPIAIJ*)A->data;
635:   Mat_SeqAIJ        *aa  =(Mat_SeqAIJ*)(mat->A)->data;
636:   Mat_SeqAIJ        *bb  =(Mat_SeqAIJ*)(mat->B)->data;

639:   ai=aa->i; aj=aa->j; adiag=aa->diag;
640:   bi=bb->i; bj=bb->j; garray = mat->garray;
641:   av=aa->a; bv=bb->a;

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

645:   if (reuse == MAT_INITIAL_MATRIX) {
646:     nza = 0;    /* num of upper triangular entries in mat->A, including diagonals */
647:     nzb = 0;    /* num of upper triangular entries in mat->B */
648:     for (i=0; i<m; i++) {
649:       nza   += (ai[i+1] - adiag[i]);
650:       countB = bi[i+1] - bi[i];
651:       bjj    = bj + bi[i];
652:       for (j=0; j<countB; j++) {
653:         if (garray[bjj[j]] > rstart) nzb++;
654:       }
655:     }

657:     nz   = nza + nzb; /* total nz of upper triangular part of mat */
658:     *nnz = nz;
659:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
660:     col  = row + nz;
661:     val  = (PetscScalar*)(col + nz);

663:     *r = row; *c = col; *v = val;
664:   } else {
665:     row = *r; col = *c; val = *v;
666:   }

668:   jj = 0; irow = rstart;
669:   for (i=0; i<m; i++) {
670:     ajj    = aj + adiag[i];                 /* ptr to the beginning of the diagonal of this row */
671:     v1     = av + adiag[i];
672:     countA = ai[i+1] - adiag[i];
673:     countB = bi[i+1] - bi[i];
674:     bjj    = bj + bi[i];
675:     v2     = bv + bi[i];

677:     /* A-part */
678:     for (j=0; j<countA; j++) {
679:       if (reuse == MAT_INITIAL_MATRIX) {
680:         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
681:       }
682:       val[jj++] = v1[j];
683:     }

685:     /* B-part */
686:     for (j=0; j < countB; j++) {
687:       if (garray[bjj[j]] > rstart) {
688:         if (reuse == MAT_INITIAL_MATRIX) {
689:           row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
690:         }
691:         val[jj++] = v2[j];
692:       }
693:     }
694:     irow++;
695:   }
696:   return(0);
697: }

699: PetscErrorCode MatDestroy_MUMPS(Mat A)
700: {
701:   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->data;

705:   PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);
706:   VecScatterDestroy(&mumps->scat_rhs);
707:   VecScatterDestroy(&mumps->scat_sol);
708:   VecDestroy(&mumps->b_seq);
709:   VecDestroy(&mumps->x_seq);
710:   PetscFree(mumps->id.perm_in);
711:   PetscFree(mumps->irn);
712:   PetscFree(mumps->info);
713:   MatMumpsResetSchur_Private(mumps);
714:   mumps->id.job = JOB_END;
715:   PetscMUMPS_c(&mumps->id);
716:   MPI_Comm_free(&mumps->comm_mumps);
717:   PetscFree(A->data);

719:   /* clear composed functions */
720:   PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverType_C",NULL);
721:   PetscObjectComposeFunction((PetscObject)A,"MatFactorSetSchurIS_C",NULL);
722:   PetscObjectComposeFunction((PetscObject)A,"MatFactorCreateSchurComplement_C",NULL);
723:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetIcntl_C",NULL);
724:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetIcntl_C",NULL);
725:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetCntl_C",NULL);
726:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetCntl_C",NULL);
727:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfo_C",NULL);
728:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfog_C",NULL);
729:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfo_C",NULL);
730:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfog_C",NULL);
731:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInverse_C",NULL);
732:   return(0);
733: }

735: PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x)
736: {
737:   Mat_MUMPS        *mumps=(Mat_MUMPS*)A->data;
738:   PetscScalar      *array;
739:   Vec              b_seq;
740:   IS               is_iden,is_petsc;
741:   PetscErrorCode   ierr;
742:   PetscInt         i;
743:   PetscBool        second_solve = PETSC_FALSE;
744:   static PetscBool cite1 = PETSC_FALSE,cite2 = PETSC_FALSE;

747:   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);
748:   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);

750:   if (A->factorerrortype) {
751:     PetscInfo2(A,"MatSolve is called with singular matrix factor, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
752:     VecSetInf(x);
753:     return(0);
754:   }

756:   mumps->id.nrhs = 1;
757:   b_seq          = mumps->b_seq;
758:   if (mumps->size > 1) {
759:     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
760:     VecScatterBegin(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);
761:     VecScatterEnd(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);
762:     if (!mumps->myid) {VecGetArray(b_seq,&array);}
763:   } else {  /* size == 1 */
764:     VecCopy(b,x);
765:     VecGetArray(x,&array);
766:   }
767:   if (!mumps->myid) { /* define rhs on the host */
768:     mumps->id.nrhs = 1;
769:     mumps->id.rhs = (MumpsScalar*)array;
770:   }

772:   /*
773:      handle condensation step of Schur complement (if any)
774:      We set by default ICNTL(26) == -1 when Schur indices have been provided by the user.
775:      According to MUMPS (5.0.0) manual, any value should be harmful during the factorization phase
776:      Unless the user provides a valid value for ICNTL(26), MatSolve and MatMatSolve routines solve the full system.
777:      This requires an extra call to PetscMUMPS_c and the computation of the factors for S
778:   */
779:   if (mumps->id.size_schur > 0 && (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2)) {
780:     if (mumps->size > 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Parallel Schur complements not yet supported from PETSc\n");
781:     second_solve = PETSC_TRUE;
782:     MatMumpsHandleSchur_Private(A,PETSC_FALSE);
783:   }
784:   /* solve phase */
785:   /*-------------*/
786:   mumps->id.job = JOB_SOLVE;
787:   PetscMUMPS_c(&mumps->id);
788:   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));

790:   /* handle expansion step of Schur complement (if any) */
791:   if (second_solve) {
792:     MatMumpsHandleSchur_Private(A,PETSC_TRUE);
793:   }

795:   if (mumps->size > 1) { /* convert mumps distributed solution to petsc mpi x */
796:     if (mumps->scat_sol && mumps->ICNTL9_pre != mumps->id.ICNTL(9)) {
797:       /* when id.ICNTL(9) changes, the contents of lsol_loc may change (not its size, lsol_loc), recreates scat_sol */
798:       VecScatterDestroy(&mumps->scat_sol);
799:     }
800:     if (!mumps->scat_sol) { /* create scatter scat_sol */
801:       ISCreateStride(PETSC_COMM_SELF,mumps->id.lsol_loc,0,1,&is_iden); /* from */
802:       for (i=0; i<mumps->id.lsol_loc; i++) {
803:         mumps->id.isol_loc[i] -= 1; /* change Fortran style to C style */
804:       }
805:       ISCreateGeneral(PETSC_COMM_SELF,mumps->id.lsol_loc,mumps->id.isol_loc,PETSC_COPY_VALUES,&is_petsc);  /* to */
806:       VecScatterCreate(mumps->x_seq,is_iden,x,is_petsc,&mumps->scat_sol);
807:       ISDestroy(&is_iden);
808:       ISDestroy(&is_petsc);

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

813:     VecScatterBegin(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);
814:     VecScatterEnd(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);
815:   }
816:   PetscLogFlops(2.0*mumps->id.RINFO(3));
817:   return(0);
818: }

820: PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x)
821: {
822:   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->data;

826:   mumps->id.ICNTL(9) = 0;
827:   MatSolve_MUMPS(A,b,x);
828:   mumps->id.ICNTL(9) = 1;
829:   return(0);
830: }

832: PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X)
833: {
835:   Mat            Bt = NULL;
836:   PetscBool      flg, flgT;
837:   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->data;
838:   PetscInt       i,nrhs,M;
839:   PetscScalar    *array,*bray;

842:   PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
843:   PetscObjectTypeCompare((PetscObject)B,MATTRANSPOSEMAT,&flgT);
844:   if (flgT) {
845:     if (mumps->size > 1) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
846:     MatTransposeGetMat(B,&Bt);
847:   } else {
848:     if (!flg) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
849:   }

851:   MatGetSize(B,&M,&nrhs);
852:   mumps->id.nrhs = nrhs;
853:   mumps->id.lrhs = M;

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

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

860:   if (mumps->size == 1) {
861:     PetscScalar *aa;
862:     PetscInt    spnr,*ia,*ja;
863:     PetscBool   second_solve = PETSC_FALSE;

865:     /* copy B to X */
866:     MatDenseGetArray(X,&array);
867:     mumps->id.rhs = (MumpsScalar*)array;
868:     if (!Bt) {
869:       MatDenseGetArray(B,&bray);
870:       PetscMemcpy(array,bray,M*nrhs*sizeof(PetscScalar));
871:       MatDenseRestoreArray(B,&bray);
872:     } else {
873:       PetscBool done;

875:       MatSeqAIJGetArray(Bt,&aa);
876:       MatGetRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&done);
877:       if (!done) SETERRQ(PetscObjectComm((PetscObject)Bt),PETSC_ERR_ARG_WRONG,"Cannot get IJ structure");
878:       mumps->id.irhs_ptr    = ia;
879:       mumps->id.irhs_sparse = ja;
880:       mumps->id.nz_rhs      = ia[spnr] - 1;
881:       mumps->id.rhs_sparse  = (MumpsScalar*)aa;
882:       mumps->id.ICNTL(20)   = 1;
883:     }
884:     /* handle condensation step of Schur complement (if any) */
885:     if (mumps->id.size_schur > 0 && (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2)) {
886:       second_solve = PETSC_TRUE;
887:       MatMumpsHandleSchur_Private(A,PETSC_FALSE);
888:     }
889:     /* solve phase */
890:     /*-------------*/
891:     mumps->id.job = JOB_SOLVE;
892:     PetscMUMPS_c(&mumps->id);
893:     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));

895:     /* handle expansion step of Schur complement (if any) */
896:     if (second_solve) {
897:       MatMumpsHandleSchur_Private(A,PETSC_TRUE);
898:     }
899:     if (Bt) {
900:       PetscBool done;

902:       MatSeqAIJRestoreArray(Bt,&aa);
903:       MatRestoreRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&done);
904:       if (!done) SETERRQ(PetscObjectComm((PetscObject)Bt),PETSC_ERR_ARG_WRONG,"Cannot restore IJ structure");
905:       mumps->id.ICNTL(20) = 0;
906:     }
907:     MatDenseRestoreArray(X,&array);
908:   } else {  /*--------- parallel case --------*/
909:     PetscInt       lsol_loc,nlsol_loc,*isol_loc,*idx,*iidx,*idxx,*isol_loc_save;
910:     MumpsScalar    *sol_loc,*sol_loc_save;
911:     IS             is_to,is_from;
912:     PetscInt       k,proc,j,m;
913:     const PetscInt *rstart;
914:     Vec            v_mpi,b_seq,x_seq;
915:     VecScatter     scat_rhs,scat_sol;

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

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

923:     lsol_loc  = mumps->id.INFO(23);
924:     nlsol_loc = nrhs*lsol_loc;     /* length of sol_loc */
925:     PetscMalloc2(nlsol_loc,&sol_loc,nlsol_loc,&isol_loc);
926:     mumps->id.sol_loc = (MumpsScalar*)sol_loc;
927:     mumps->id.isol_loc = isol_loc;

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

931:     /* copy rhs matrix B into vector v_mpi */
932:     MatGetLocalSize(B,&m,NULL);
933:     MatDenseGetArray(B,&bray);
934:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)B),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi);
935:     MatDenseRestoreArray(B,&bray);

937:     /* scatter v_mpi to b_seq because MUMPS only supports centralized rhs */
938:     /* idx: maps from k-th index of v_mpi to (i,j)-th global entry of B;
939:       iidx: inverse of idx, will be used by scattering xx_seq -> X       */
940:     PetscMalloc2(nrhs*M,&idx,nrhs*M,&iidx);
941:     MatGetOwnershipRanges(B,&rstart);
942:     k = 0;
943:     for (proc=0; proc<mumps->size; proc++){
944:       for (j=0; j<nrhs; j++){
945:         for (i=rstart[proc]; i<rstart[proc+1]; i++){
946:           iidx[j*M + i] = k;
947:           idx[k++]      = j*M + i;
948:         }
949:       }
950:     }

952:     if (!mumps->myid) {
953:       VecCreateSeq(PETSC_COMM_SELF,nrhs*M,&b_seq);
954:       ISCreateGeneral(PETSC_COMM_SELF,nrhs*M,idx,PETSC_COPY_VALUES,&is_to);
955:       ISCreateStride(PETSC_COMM_SELF,nrhs*M,0,1,&is_from);
956:     } else {
957:       VecCreateSeq(PETSC_COMM_SELF,0,&b_seq);
958:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_to);
959:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_from);
960:     }
961:     VecScatterCreate(v_mpi,is_from,b_seq,is_to,&scat_rhs);
962:     VecScatterBegin(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);
963:     ISDestroy(&is_to);
964:     ISDestroy(&is_from);
965:     VecScatterEnd(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);

967:     if (!mumps->myid) { /* define rhs on the host */
968:       VecGetArray(b_seq,&bray);
969:       mumps->id.rhs = (MumpsScalar*)bray;
970:       VecRestoreArray(b_seq,&bray);
971:     }

973:     /* solve phase */
974:     /*-------------*/
975:     mumps->id.job = JOB_SOLVE;
976:     PetscMUMPS_c(&mumps->id);
977:     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));

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

983:     /* create scatter scat_sol */
984:     PetscMalloc1(nlsol_loc,&idxx);
985:     ISCreateStride(PETSC_COMM_SELF,nlsol_loc,0,1,&is_from);
986:     for (i=0; i<lsol_loc; i++) {
987:       isol_loc[i] -= 1; /* change Fortran style to C style */
988:       idxx[i] = iidx[isol_loc[i]];
989:       for (j=1; j<nrhs; j++){
990:         idxx[j*lsol_loc+i] = iidx[isol_loc[i]+j*M];
991:       }
992:     }
993:     ISCreateGeneral(PETSC_COMM_SELF,nlsol_loc,idxx,PETSC_COPY_VALUES,&is_to);
994:     VecScatterCreate(x_seq,is_from,v_mpi,is_to,&scat_sol);
995:     VecScatterBegin(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);
996:     ISDestroy(&is_from);
997:     ISDestroy(&is_to);
998:     VecScatterEnd(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);
999:     MatDenseRestoreArray(X,&array);

1001:     /* free spaces */
1002:     mumps->id.sol_loc = sol_loc_save;
1003:     mumps->id.isol_loc = isol_loc_save;

1005:     PetscFree2(sol_loc,isol_loc);
1006:     PetscFree2(idx,iidx);
1007:     PetscFree(idxx);
1008:     VecDestroy(&x_seq);
1009:     VecDestroy(&v_mpi);
1010:     VecDestroy(&b_seq);
1011:     VecScatterDestroy(&scat_rhs);
1012:     VecScatterDestroy(&scat_sol);
1013:   }
1014:   PetscLogFlops(2.0*nrhs*mumps->id.RINFO(3));
1015:   return(0);
1016: }

1018: #if !defined(PETSC_USE_COMPLEX)
1019: /*
1020:   input:
1021:    F:        numeric factor
1022:   output:
1023:    nneg:     total number of negative pivots
1024:    nzero:    total number of zero pivots
1025:    npos:     (global dimension of F) - nneg - nzero
1026: */
1027: PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
1028: {
1029:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1031:   PetscMPIInt    size;

1034:   MPI_Comm_size(PetscObjectComm((PetscObject)F),&size);
1035:   /* 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 */
1036:   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));

1038:   if (nneg) *nneg = mumps->id.INFOG(12);
1039:   if (nzero || npos) {
1040:     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");
1041:     if (nzero) *nzero = mumps->id.INFOG(28);
1042:     if (npos) *npos   = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28));
1043:   }
1044:   return(0);
1045: }
1046: #endif

1048: PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
1049: {
1050:   Mat_MUMPS      *mumps =(Mat_MUMPS*)(F)->data;
1052:   PetscBool      isMPIAIJ;

1055:   if (mumps->id.INFOG(1) < 0) {
1056:     if (mumps->id.INFOG(1) == -6) {
1057:       PetscInfo2(A,"MatFactorNumeric is called with singular matrix structure, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1058:     }
1059:     PetscInfo2(A,"MatFactorNumeric is called after analysis phase fails, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1060:     return(0);
1061:   }

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

1065:   /* numerical factorization phase */
1066:   /*-------------------------------*/
1067:   mumps->id.job = JOB_FACTNUMERIC;
1068:   if (!mumps->id.ICNTL(18)) { /* A is centralized */
1069:     if (!mumps->myid) {
1070:       mumps->id.a = (MumpsScalar*)mumps->val;
1071:     }
1072:   } else {
1073:     mumps->id.a_loc = (MumpsScalar*)mumps->val;
1074:   }
1075:   PetscMUMPS_c(&mumps->id);
1076:   if (mumps->id.INFOG(1) < 0) {
1077:     if (A->erroriffailure) {
1078:       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));
1079:     } else {
1080:       if (mumps->id.INFOG(1) == -10) { /* numerically singular matrix */
1081:         PetscInfo2(F,"matrix is numerically singular, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1082:         F->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1083:       } else if (mumps->id.INFOG(1) == -13) {
1084:         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));
1085:         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1086:       } else if (mumps->id.INFOG(1) == -8 || mumps->id.INFOG(1) == -9 || (-16 < mumps->id.INFOG(1) && mumps->id.INFOG(1) < -10) ) {
1087:         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));
1088:         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1089:       } else {
1090:         PetscInfo2(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1091:         F->factorerrortype = MAT_FACTOR_OTHER;
1092:       }
1093:     }
1094:   }
1095:   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));

1097:   F->assembled    = PETSC_TRUE;
1098:   mumps->matstruc = SAME_NONZERO_PATTERN;
1099:   if (F->schur) { /* reset Schur status to unfactored */
1100:     if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
1101:       mumps->id.ICNTL(19) = 2;
1102:       MatTranspose(F->schur,MAT_INPLACE_MATRIX,&F->schur);
1103:     }
1104:     MatFactorRestoreSchurComplement(F,NULL,MAT_FACTOR_SCHUR_UNFACTORED);
1105:   }

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

1110:   if (mumps->size > 1) {
1111:     PetscInt    lsol_loc;
1112:     PetscScalar *sol_loc;

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

1116:     /* distributed solution; Create x_seq=sol_loc for repeated use */
1117:     if (mumps->x_seq) {
1118:       VecScatterDestroy(&mumps->scat_sol);
1119:       PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);
1120:       VecDestroy(&mumps->x_seq);
1121:     }
1122:     lsol_loc = mumps->id.INFO(23); /* length of sol_loc */
1123:     PetscMalloc2(lsol_loc,&sol_loc,lsol_loc,&mumps->id.isol_loc);
1124:     mumps->id.lsol_loc = lsol_loc;
1125:     mumps->id.sol_loc = (MumpsScalar*)sol_loc;
1126:     VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq);
1127:   }
1128:   PetscLogFlops(mumps->id.RINFO(2));
1129:   return(0);
1130: }

1132: /* Sets MUMPS options from the options database */
1133: PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A)
1134: {
1135:   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1137:   PetscInt       icntl,info[40],i,ninfo=40;
1138:   PetscBool      flg;

1141:   PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MUMPS Options","Mat");
1142:   PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);
1143:   if (flg) mumps->id.ICNTL(1) = icntl;
1144:   PetscOptionsInt("-mat_mumps_icntl_2","ICNTL(2): output stream for diagnostic printing, statistics, and warning","None",mumps->id.ICNTL(2),&icntl,&flg);
1145:   if (flg) mumps->id.ICNTL(2) = icntl;
1146:   PetscOptionsInt("-mat_mumps_icntl_3","ICNTL(3): output stream for global information, collected on the host","None",mumps->id.ICNTL(3),&icntl,&flg);
1147:   if (flg) mumps->id.ICNTL(3) = icntl;

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

1153:   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);
1154:   if (flg) mumps->id.ICNTL(6) = icntl;

1156:   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);
1157:   if (flg) {
1158:     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");
1159:     else mumps->id.ICNTL(7) = icntl;
1160:   }

1162:   PetscOptionsInt("-mat_mumps_icntl_8","ICNTL(8): scaling strategy (-2 to 8 or 77)","None",mumps->id.ICNTL(8),&mumps->id.ICNTL(8),NULL);
1163:   /* 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() */
1164:   PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",mumps->id.ICNTL(10),&mumps->id.ICNTL(10),NULL);
1165:   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);
1166:   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);
1167:   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);
1168:   PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage increase in the estimated working space","None",mumps->id.ICNTL(14),&mumps->id.ICNTL(14),NULL);
1169:   PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): computes the Schur complement","None",mumps->id.ICNTL(19),&mumps->id.ICNTL(19),NULL);
1170:   if (mumps->id.ICNTL(19) <= 0 || mumps->id.ICNTL(19) > 3) { /* reset any schur data (if any) */
1171:     MatDestroy(&F->schur);
1172:     MatMumpsResetSchur_Private(mumps);
1173:   }
1174:   /* 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 */
1175:   /* 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 */

1177:   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);
1178:   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);
1179:   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);
1180:   if (mumps->id.ICNTL(24)) {
1181:     mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */
1182:   }

1184:   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);
1185:   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);
1186:   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);
1187:   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);
1188:   PetscOptionsInt("-mat_mumps_icntl_29","ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis","None",mumps->id.ICNTL(29),&mumps->id.ICNTL(29),NULL);
1189:   /* 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 */
1190:   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);
1191:   /* 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 */
1192:   PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),NULL);
1193:   PetscOptionsInt("-mat_mumps_icntl_35","ICNTL(35): activates Block Lock Rank (BLR) based factorization","None",mumps->id.ICNTL(35),&mumps->id.ICNTL(35),NULL);

1195:   PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),NULL);
1196:   PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",mumps->id.CNTL(2),&mumps->id.CNTL(2),NULL);
1197:   PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),NULL);
1198:   PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",mumps->id.CNTL(4),&mumps->id.CNTL(4),NULL);
1199:   PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",mumps->id.CNTL(5),&mumps->id.CNTL(5),NULL);
1200:   PetscOptionsReal("-mat_mumps_cntl_7","CNTL(7): dropping parameter used during BLR","None",mumps->id.CNTL(7),&mumps->id.CNTL(7),NULL);

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

1204:   PetscOptionsIntArray("-mat_mumps_view_info","request INFO local to each processor","",info,&ninfo,NULL);
1205:   if (ninfo) {
1206:     if (ninfo > 40) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"number of INFO %d must <= 40\n",ninfo);
1207:     PetscMalloc1(ninfo,&mumps->info);
1208:     mumps->ninfo = ninfo;
1209:     for (i=0; i<ninfo; i++) {
1210:       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);
1211:       else  mumps->info[i] = info[i];
1212:     }
1213:   }

1215:   PetscOptionsEnd();
1216:   return(0);
1217: }

1219: PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps)
1220: {

1224:   MPI_Comm_rank(PetscObjectComm((PetscObject)A), &mumps->myid);
1225:   MPI_Comm_size(PetscObjectComm((PetscObject)A),&mumps->size);
1226:   MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(mumps->comm_mumps));

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

1230:   mumps->id.job = JOB_INIT;
1231:   mumps->id.par = 1;  /* host participates factorizaton and solve */
1232:   mumps->id.sym = mumps->sym;
1233:   PetscMUMPS_c(&mumps->id);

1235:   mumps->scat_rhs     = NULL;
1236:   mumps->scat_sol     = NULL;

1238:   /* set PETSc-MUMPS default options - override MUMPS default */
1239:   mumps->id.ICNTL(3) = 0;
1240:   mumps->id.ICNTL(4) = 0;
1241:   if (mumps->size == 1) {
1242:     mumps->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
1243:   } else {
1244:     mumps->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
1245:     mumps->id.ICNTL(20) = 0;   /* rhs is in dense format */
1246:     mumps->id.ICNTL(21) = 1;   /* distributed solution */
1247:   }

1249:   /* schur */
1250:   mumps->id.size_schur      = 0;
1251:   mumps->id.listvar_schur   = NULL;
1252:   mumps->id.schur           = NULL;
1253:   mumps->sizeredrhs         = 0;
1254:   mumps->schur_sol          = NULL;
1255:   mumps->schur_sizesol      = 0;
1256:   return(0);
1257: }

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

1264:   if (mumps->id.INFOG(1) < 0) {
1265:     if (A->erroriffailure) {
1266:       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
1267:     } else {
1268:       if (mumps->id.INFOG(1) == -6) {
1269:         PetscInfo2(F,"matrix is singular in structure, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1270:         F->factorerrortype = MAT_FACTOR_STRUCT_ZEROPIVOT;
1271:       } else if (mumps->id.INFOG(1) == -5 || mumps->id.INFOG(1) == -7) {
1272:         PetscInfo2(F,"problem of workspace, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1273:         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1274:       } else {
1275:         PetscInfo2(F,"Error reported by MUMPS in analysis phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1276:         F->factorerrortype = MAT_FACTOR_OTHER;
1277:       }
1278:     }
1279:   }
1280:   return(0);
1281: }

1283: /* Note Petsc r(=c) permutation is used when mumps->id.ICNTL(7)==1 with centralized assembled matrix input; otherwise r and c are ignored */
1284: PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1285: {
1286:   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1288:   Vec            b;
1289:   IS             is_iden;
1290:   const PetscInt M = A->rmap->N;

1293:   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;

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

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

1300:   /* analysis phase */
1301:   /*----------------*/
1302:   mumps->id.job = JOB_FACTSYMBOLIC;
1303:   mumps->id.n   = M;
1304:   switch (mumps->id.ICNTL(18)) {
1305:   case 0:  /* centralized assembled matrix input */
1306:     if (!mumps->myid) {
1307:       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1308:       if (mumps->id.ICNTL(6)>1) {
1309:         mumps->id.a = (MumpsScalar*)mumps->val;
1310:       }
1311:       if (mumps->id.ICNTL(7) == 1) { /* use user-provide matrix ordering - assuming r = c ordering */
1312:         /*
1313:         PetscBool      flag;
1314:         ISEqual(r,c,&flag);
1315:         if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"row_perm != col_perm");
1316:         ISView(r,PETSC_VIEWER_STDOUT_SELF);
1317:          */
1318:         if (!mumps->myid) {
1319:           const PetscInt *idx;
1320:           PetscInt       i,*perm_in;

1322:           PetscMalloc1(M,&perm_in);
1323:           ISGetIndices(r,&idx);

1325:           mumps->id.perm_in = perm_in;
1326:           for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */
1327:           ISRestoreIndices(r,&idx);
1328:         }
1329:       }
1330:     }
1331:     break;
1332:   case 3:  /* distributed assembled matrix input (size>1) */
1333:     mumps->id.nz_loc = mumps->nz;
1334:     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1335:     if (mumps->id.ICNTL(6)>1) {
1336:       mumps->id.a_loc = (MumpsScalar*)mumps->val;
1337:     }
1338:     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1339:     if (!mumps->myid) {
1340:       VecCreateSeq(PETSC_COMM_SELF,A->rmap->N,&mumps->b_seq);
1341:       ISCreateStride(PETSC_COMM_SELF,A->rmap->N,0,1,&is_iden);
1342:     } else {
1343:       VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);
1344:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);
1345:     }
1346:     MatCreateVecs(A,NULL,&b);
1347:     VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);
1348:     ISDestroy(&is_iden);
1349:     VecDestroy(&b);
1350:     break;
1351:   }
1352:   PetscMUMPS_c(&mumps->id);
1353:   MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);

1355:   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1356:   F->ops->solve           = MatSolve_MUMPS;
1357:   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1358:   F->ops->matsolve        = MatMatSolve_MUMPS;
1359:   return(0);
1360: }

1362: /* Note the Petsc r and c permutations are ignored */
1363: PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1364: {
1365:   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1367:   Vec            b;
1368:   IS             is_iden;
1369:   const PetscInt M = A->rmap->N;

1372:   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;

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

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

1379:   /* analysis phase */
1380:   /*----------------*/
1381:   mumps->id.job = JOB_FACTSYMBOLIC;
1382:   mumps->id.n   = M;
1383:   switch (mumps->id.ICNTL(18)) {
1384:   case 0:  /* centralized assembled matrix input */
1385:     if (!mumps->myid) {
1386:       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1387:       if (mumps->id.ICNTL(6)>1) {
1388:         mumps->id.a = (MumpsScalar*)mumps->val;
1389:       }
1390:     }
1391:     break;
1392:   case 3:  /* distributed assembled matrix input (size>1) */
1393:     mumps->id.nz_loc = mumps->nz;
1394:     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1395:     if (mumps->id.ICNTL(6)>1) {
1396:       mumps->id.a_loc = (MumpsScalar*)mumps->val;
1397:     }
1398:     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1399:     if (!mumps->myid) {
1400:       VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);
1401:       ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);
1402:     } else {
1403:       VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);
1404:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);
1405:     }
1406:     MatCreateVecs(A,NULL,&b);
1407:     VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);
1408:     ISDestroy(&is_iden);
1409:     VecDestroy(&b);
1410:     break;
1411:   }
1412:   PetscMUMPS_c(&mumps->id);
1413:   MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);

1415:   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1416:   F->ops->solve           = MatSolve_MUMPS;
1417:   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1418:   return(0);
1419: }

1421: /* Note the Petsc r permutation and factor info are ignored */
1422: PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
1423: {
1424:   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1426:   Vec            b;
1427:   IS             is_iden;
1428:   const PetscInt M = A->rmap->N;

1431:   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;

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

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

1438:   /* analysis phase */
1439:   /*----------------*/
1440:   mumps->id.job = JOB_FACTSYMBOLIC;
1441:   mumps->id.n   = M;
1442:   switch (mumps->id.ICNTL(18)) {
1443:   case 0:  /* centralized assembled matrix input */
1444:     if (!mumps->myid) {
1445:       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1446:       if (mumps->id.ICNTL(6)>1) {
1447:         mumps->id.a = (MumpsScalar*)mumps->val;
1448:       }
1449:     }
1450:     break;
1451:   case 3:  /* distributed assembled matrix input (size>1) */
1452:     mumps->id.nz_loc = mumps->nz;
1453:     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1454:     if (mumps->id.ICNTL(6)>1) {
1455:       mumps->id.a_loc = (MumpsScalar*)mumps->val;
1456:     }
1457:     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1458:     if (!mumps->myid) {
1459:       VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);
1460:       ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);
1461:     } else {
1462:       VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);
1463:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);
1464:     }
1465:     MatCreateVecs(A,NULL,&b);
1466:     VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);
1467:     ISDestroy(&is_iden);
1468:     VecDestroy(&b);
1469:     break;
1470:   }
1471:   PetscMUMPS_c(&mumps->id);
1472:   MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);

1474:   F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
1475:   F->ops->solve                 = MatSolve_MUMPS;
1476:   F->ops->solvetranspose        = MatSolve_MUMPS;
1477:   F->ops->matsolve              = MatMatSolve_MUMPS;
1478: #if defined(PETSC_USE_COMPLEX)
1479:   F->ops->getinertia = NULL;
1480: #else
1481:   F->ops->getinertia = MatGetInertia_SBAIJMUMPS;
1482: #endif
1483:   return(0);
1484: }

1486: PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
1487: {
1488:   PetscErrorCode    ierr;
1489:   PetscBool         iascii;
1490:   PetscViewerFormat format;
1491:   Mat_MUMPS         *mumps=(Mat_MUMPS*)A->data;

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

1497:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1498:   if (iascii) {
1499:     PetscViewerGetFormat(viewer,&format);
1500:     if (format == PETSC_VIEWER_ASCII_INFO) {
1501:       PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");
1502:       PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                   %d \n",mumps->id.sym);
1503:       PetscViewerASCIIPrintf(viewer,"  PAR (host participation):            %d \n",mumps->id.par);
1504:       PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):         %d \n",mumps->id.ICNTL(1));
1505:       PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg): %d \n",mumps->id.ICNTL(2));
1506:       PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):   %d \n",mumps->id.ICNTL(3));
1507:       PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):        %d \n",mumps->id.ICNTL(4));
1508:       PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):         %d \n",mumps->id.ICNTL(5));
1509:       PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):        %d \n",mumps->id.ICNTL(6));
1510:       PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (sequential matrix ordering):%d \n",mumps->id.ICNTL(7));
1511:       PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scaling strategy):        %d \n",mumps->id.ICNTL(8));
1512:       PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements):  %d \n",mumps->id.ICNTL(10));
1513:       PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):          %d \n",mumps->id.ICNTL(11));
1514:       if (mumps->id.ICNTL(11)>0) {
1515:         PetscViewerASCIIPrintf(viewer,"    RINFOG(4) (inf norm of input mat):        %g\n",mumps->id.RINFOG(4));
1516:         PetscViewerASCIIPrintf(viewer,"    RINFOG(5) (inf norm of solution):         %g\n",mumps->id.RINFOG(5));
1517:         PetscViewerASCIIPrintf(viewer,"    RINFOG(6) (inf norm of residual):         %g\n",mumps->id.RINFOG(6));
1518:         PetscViewerASCIIPrintf(viewer,"    RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8));
1519:         PetscViewerASCIIPrintf(viewer,"    RINFOG(9) (error estimate):               %g \n",mumps->id.RINFOG(9));
1520:         PetscViewerASCIIPrintf(viewer,"    RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11));
1521:       }
1522:       PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",mumps->id.ICNTL(12));
1523:       PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",mumps->id.ICNTL(13));
1524:       PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",mumps->id.ICNTL(14));
1525:       /* ICNTL(15-17) not used */
1526:       PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",mumps->id.ICNTL(18));
1527:       PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Schur complement info):                      %d \n",mumps->id.ICNTL(19));
1528:       PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",mumps->id.ICNTL(20));
1529:       PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d \n",mumps->id.ICNTL(21));
1530:       PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",mumps->id.ICNTL(22));
1531:       PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",mumps->id.ICNTL(23));

1533:       PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",mumps->id.ICNTL(24));
1534:       PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",mumps->id.ICNTL(25));
1535:       PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",mumps->id.ICNTL(26));
1536:       PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",mumps->id.ICNTL(27));
1537:       PetscViewerASCIIPrintf(viewer,"  ICNTL(28) (use parallel or sequential ordering):        %d \n",mumps->id.ICNTL(28));
1538:       PetscViewerASCIIPrintf(viewer,"  ICNTL(29) (parallel ordering):                          %d \n",mumps->id.ICNTL(29));

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

1545:       PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",mumps->id.CNTL(1));
1546:       PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2));
1547:       PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",mumps->id.CNTL(3));
1548:       PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",mumps->id.CNTL(4));
1549:       PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",mumps->id.CNTL(5));
1550:       PetscViewerASCIIPrintf(viewer,"  CNTL(7) (dropping parameter for BLR):       %g \n",mumps->id.CNTL(7));

1552:       /* infomation local to each processor */
1553:       PetscViewerASCIIPrintf(viewer, "  RINFO(1) (local estimated flops for the elimination after analysis): \n");
1554:       PetscViewerASCIIPushSynchronized(viewer);
1555:       PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %g \n",mumps->myid,mumps->id.RINFO(1));
1556:       PetscViewerFlush(viewer);
1557:       PetscViewerASCIIPrintf(viewer, "  RINFO(2) (local estimated flops for the assembly after factorization): \n");
1558:       PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(2));
1559:       PetscViewerFlush(viewer);
1560:       PetscViewerASCIIPrintf(viewer, "  RINFO(3) (local estimated flops for the elimination after factorization): \n");
1561:       PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(3));
1562:       PetscViewerFlush(viewer);

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

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

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

1576:       if (mumps->ninfo && mumps->ninfo <= 40){
1577:         PetscInt i;
1578:         for (i=0; i<mumps->ninfo; i++){
1579:           PetscViewerASCIIPrintf(viewer, "  INFO(%d): \n",mumps->info[i]);
1580:           PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(mumps->info[i]));
1581:           PetscViewerFlush(viewer);
1582:         }
1583:       }
1584:       PetscViewerASCIIPopSynchronized(viewer);

1586:       if (!mumps->myid) { /* information from the host */
1587:         PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",mumps->id.RINFOG(1));
1588:         PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",mumps->id.RINFOG(2));
1589:         PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",mumps->id.RINFOG(3));
1590:         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));

1592:         PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(3));
1593:         PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(4));
1594:         PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",mumps->id.INFOG(5));
1595:         PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",mumps->id.INFOG(6));
1596:         PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively use after analysis): %d \n",mumps->id.INFOG(7));
1597:         PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",mumps->id.INFOG(8));
1598:         PetscViewerASCIIPrintf(viewer,"  INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",mumps->id.INFOG(9));
1599:         PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",mumps->id.INFOG(10));
1600:         PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",mumps->id.INFOG(11));
1601:         PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",mumps->id.INFOG(12));
1602:         PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",mumps->id.INFOG(13));
1603:         PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",mumps->id.INFOG(14));
1604:         PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",mumps->id.INFOG(15));
1605:         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));
1606:         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));
1607:         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));
1608:         PetscViewerASCIIPrintf(viewer,"  INFOG(19) (size of all MUMPS internal data allocated during factorization: sum over all processors): %d \n",mumps->id.INFOG(19));
1609:         PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",mumps->id.INFOG(20));
1610:         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));
1611:         PetscViewerASCIIPrintf(viewer,"  INFOG(22) (size in MB of memory effectively used during factorization - sum over all processors): %d \n",mumps->id.INFOG(22));
1612:         PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",mumps->id.INFOG(23));
1613:         PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",mumps->id.INFOG(24));
1614:         PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",mumps->id.INFOG(25));
1615:         PetscViewerASCIIPrintf(viewer,"  INFOG(28) (after factorization: number of null pivots encountered): %d\n",mumps->id.INFOG(28));
1616:         PetscViewerASCIIPrintf(viewer,"  INFOG(29) (after factorization: effective number of entries in the factors (sum over all processors)): %d\n",mumps->id.INFOG(29));
1617:         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));
1618:         PetscViewerASCIIPrintf(viewer,"  INFOG(32) (after analysis: type of analysis done): %d\n",mumps->id.INFOG(32));
1619:         PetscViewerASCIIPrintf(viewer,"  INFOG(33) (value used for ICNTL(8)): %d\n",mumps->id.INFOG(33));
1620:         PetscViewerASCIIPrintf(viewer,"  INFOG(34) (exponent of the determinant if determinant is requested): %d\n",mumps->id.INFOG(34));
1621:       }
1622:     }
1623:   }
1624:   return(0);
1625: }

1627: PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
1628: {
1629:   Mat_MUMPS *mumps =(Mat_MUMPS*)A->data;

1632:   info->block_size        = 1.0;
1633:   info->nz_allocated      = mumps->id.INFOG(20);
1634:   info->nz_used           = mumps->id.INFOG(20);
1635:   info->nz_unneeded       = 0.0;
1636:   info->assemblies        = 0.0;
1637:   info->mallocs           = 0.0;
1638:   info->memory            = 0.0;
1639:   info->fill_ratio_given  = 0;
1640:   info->fill_ratio_needed = 0;
1641:   info->factor_mallocs    = 0;
1642:   return(0);
1643: }

1645: /* -------------------------------------------------------------------------------------------*/
1646: PetscErrorCode MatFactorSetSchurIS_MUMPS(Mat F, IS is)
1647: {
1648:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1649:   const PetscInt *idxs;
1650:   PetscInt       size,i;

1654:   ISGetLocalSize(is,&size);
1655:   if (mumps->size > 1) {
1656:     PetscBool ls,gs;

1658:     ls   = mumps->myid ? (size ? PETSC_FALSE : PETSC_TRUE) : PETSC_TRUE;
1659:     MPI_Allreduce(&ls,&gs,1,MPIU_BOOL,MPI_LAND,mumps->comm_mumps);
1660:     if (!gs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MUMPS distributed parallel Schur complements not yet supported from PETSc\n");
1661:   }
1662:   if (mumps->id.size_schur != size) {
1663:     PetscFree2(mumps->id.listvar_schur,mumps->id.schur);
1664:     mumps->id.size_schur = size;
1665:     mumps->id.schur_lld  = size;
1666:     PetscMalloc2(size,&mumps->id.listvar_schur,size*size,&mumps->id.schur);
1667:   }

1669:   /* Schur complement matrix */
1670:   MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.size_schur,(PetscScalar*)mumps->id.schur,&F->schur);
1671:   if (mumps->sym == 1) {
1672:     MatSetOption(F->schur,MAT_SPD,PETSC_TRUE);
1673:   }

1675:   /* MUMPS expects Fortran style indices */
1676:   ISGetIndices(is,&idxs);
1677:   PetscMemcpy(mumps->id.listvar_schur,idxs,size*sizeof(PetscInt));
1678:   for (i=0;i<size;i++) mumps->id.listvar_schur[i]++;
1679:   ISRestoreIndices(is,&idxs);
1680:   if (mumps->size > 1) {
1681:     mumps->id.ICNTL(19) = 1; /* MUMPS returns Schur centralized on the host */
1682:   } else {
1683:     if (F->factortype == MAT_FACTOR_LU) {
1684:       mumps->id.ICNTL(19) = 3; /* MUMPS returns full matrix */
1685:     } else {
1686:       mumps->id.ICNTL(19) = 2; /* MUMPS returns lower triangular part */
1687:     }
1688:   }
1689:   /* set a special value of ICNTL (not handled my MUMPS) to be used in the solve phase by PETSc */
1690:   mumps->id.ICNTL(26) = -1;
1691:   return(0);
1692: }

1694: /* -------------------------------------------------------------------------------------------*/
1695: PetscErrorCode MatFactorCreateSchurComplement_MUMPS(Mat F,Mat* S)
1696: {
1697:   Mat            St;
1698:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1699:   PetscScalar    *array;
1700: #if defined(PETSC_USE_COMPLEX)
1701:   PetscScalar    im = PetscSqrtScalar((PetscScalar)-1.0);
1702: #endif

1706:   if (!mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it");
1707:   MatCreate(PETSC_COMM_SELF,&St);
1708:   MatSetSizes(St,PETSC_DECIDE,PETSC_DECIDE,mumps->id.size_schur,mumps->id.size_schur);
1709:   MatSetType(St,MATDENSE);
1710:   MatSetUp(St);
1711:   MatDenseGetArray(St,&array);
1712:   if (!mumps->sym) { /* MUMPS always return a full matrix */
1713:     if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
1714:       PetscInt i,j,N=mumps->id.size_schur;
1715:       for (i=0;i<N;i++) {
1716:         for (j=0;j<N;j++) {
1717: #if !defined(PETSC_USE_COMPLEX)
1718:           PetscScalar val = mumps->id.schur[i*N+j];
1719: #else
1720:           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
1721: #endif
1722:           array[j*N+i] = val;
1723:         }
1724:       }
1725:     } else { /* stored by columns */
1726:       PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));
1727:     }
1728:   } else { /* either full or lower-triangular (not packed) */
1729:     if (mumps->id.ICNTL(19) == 2) { /* lower triangular stored by columns */
1730:       PetscInt i,j,N=mumps->id.size_schur;
1731:       for (i=0;i<N;i++) {
1732:         for (j=i;j<N;j++) {
1733: #if !defined(PETSC_USE_COMPLEX)
1734:           PetscScalar val = mumps->id.schur[i*N+j];
1735: #else
1736:           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
1737: #endif
1738:           array[i*N+j] = val;
1739:           array[j*N+i] = val;
1740:         }
1741:       }
1742:     } else if (mumps->id.ICNTL(19) == 3) { /* full matrix */
1743:       PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));
1744:     } else { /* ICNTL(19) == 1 lower triangular stored by rows */
1745:       PetscInt i,j,N=mumps->id.size_schur;
1746:       for (i=0;i<N;i++) {
1747:         for (j=0;j<i+1;j++) {
1748: #if !defined(PETSC_USE_COMPLEX)
1749:           PetscScalar val = mumps->id.schur[i*N+j];
1750: #else
1751:           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
1752: #endif
1753:           array[i*N+j] = val;
1754:           array[j*N+i] = val;
1755:         }
1756:       }
1757:     }
1758:   }
1759:   MatDenseRestoreArray(St,&array);
1760:   *S   = St;
1761:   return(0);
1762: }

1764: /* -------------------------------------------------------------------------------------------*/
1765: PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival)
1766: {
1767:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

1770:   mumps->id.ICNTL(icntl) = ival;
1771:   return(0);
1772: }

1774: PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt *ival)
1775: {
1776:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

1779:   *ival = mumps->id.ICNTL(icntl);
1780:   return(0);
1781: }

1783: /*@
1784:   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()

1786:    Logically Collective on Mat

1788:    Input Parameters:
1789: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1790: .  icntl - index of MUMPS parameter array ICNTL()
1791: -  ival - value of MUMPS ICNTL(icntl)

1793:   Options Database:
1794: .   -mat_mumps_icntl_<icntl> <ival>

1796:    Level: beginner

1798:    References:
1799: .     MUMPS Users' Guide

1801: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
1802:  @*/
1803: PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival)
1804: {
1806: 
1809:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
1812:   PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));
1813:   return(0);
1814: }

1816: /*@
1817:   MatMumpsGetIcntl - Get MUMPS parameter ICNTL()

1819:    Logically Collective on Mat

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

1825:   Output Parameter:
1826: .  ival - value of MUMPS ICNTL(icntl)

1828:    Level: beginner

1830:    References:
1831: .     MUMPS Users' Guide

1833: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
1834: @*/
1835: PetscErrorCode MatMumpsGetIcntl(Mat F,PetscInt icntl,PetscInt *ival)
1836: {

1841:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
1844:   PetscUseMethod(F,"MatMumpsGetIcntl_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));
1845:   return(0);
1846: }

1848: /* -------------------------------------------------------------------------------------------*/
1849: PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal val)
1850: {
1851:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

1854:   mumps->id.CNTL(icntl) = val;
1855:   return(0);
1856: }

1858: PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal *val)
1859: {
1860:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

1863:   *val = mumps->id.CNTL(icntl);
1864:   return(0);
1865: }

1867: /*@
1868:   MatMumpsSetCntl - Set MUMPS parameter CNTL()

1870:    Logically Collective on Mat

1872:    Input Parameters:
1873: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1874: .  icntl - index of MUMPS parameter array CNTL()
1875: -  val - value of MUMPS CNTL(icntl)

1877:   Options Database:
1878: .   -mat_mumps_cntl_<icntl> <val>

1880:    Level: beginner

1882:    References:
1883: .     MUMPS Users' Guide

1885: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
1886: @*/
1887: PetscErrorCode MatMumpsSetCntl(Mat F,PetscInt icntl,PetscReal val)
1888: {

1893:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
1896:   PetscTryMethod(F,"MatMumpsSetCntl_C",(Mat,PetscInt,PetscReal),(F,icntl,val));
1897:   return(0);
1898: }

1900: /*@
1901:   MatMumpsGetCntl - Get MUMPS parameter CNTL()

1903:    Logically Collective on Mat

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

1909:   Output Parameter:
1910: .  val - value of MUMPS CNTL(icntl)

1912:    Level: beginner

1914:    References:
1915: .      MUMPS Users' Guide

1917: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
1918: @*/
1919: PetscErrorCode MatMumpsGetCntl(Mat F,PetscInt icntl,PetscReal *val)
1920: {

1925:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
1928:   PetscUseMethod(F,"MatMumpsGetCntl_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));
1929:   return(0);
1930: }

1932: PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F,PetscInt icntl,PetscInt *info)
1933: {
1934:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

1937:   *info = mumps->id.INFO(icntl);
1938:   return(0);
1939: }

1941: PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F,PetscInt icntl,PetscInt *infog)
1942: {
1943:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

1946:   *infog = mumps->id.INFOG(icntl);
1947:   return(0);
1948: }

1950: PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfo)
1951: {
1952:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

1955:   *rinfo = mumps->id.RINFO(icntl);
1956:   return(0);
1957: }

1959: PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfog)
1960: {
1961:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

1964:   *rinfog = mumps->id.RINFOG(icntl);
1965:   return(0);
1966: }

1968: PetscErrorCode MatMumpsGetInverse_MUMPS(Mat F,Mat spRHS)
1969: {
1971:   Mat            Bt = NULL;
1972:   PetscBool      flgT;
1973:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1974:   PetscBool      done;
1975:   PetscScalar    *aa;
1976:   PetscInt       spnr,*ia,*ja;

1979:   if (!mumps->myid) {
1981:     PetscObjectTypeCompare((PetscObject)spRHS,MATTRANSPOSEMAT,&flgT);
1982:     if (flgT) {
1983:       MatTransposeGetMat(spRHS,&Bt);
1984:     } else {
1985:       SETERRQ(PetscObjectComm((PetscObject)spRHS),PETSC_ERR_ARG_WRONG,"Matrix spRHS must be type MATTRANSPOSEMAT matrix");
1986:     }
1987:   }

1989:   MatMumpsSetIcntl(F,30,1);

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

1996:     mumps->id.irhs_ptr    = ia;
1997:     mumps->id.irhs_sparse = ja;
1998:     mumps->id.nz_rhs      = ia[spnr] - 1;
1999:     mumps->id.rhs_sparse  = (MumpsScalar*)aa;
2000:   } else {
2001:     mumps->id.irhs_ptr    = NULL;
2002:     mumps->id.irhs_sparse = NULL;
2003:     mumps->id.nz_rhs      = 0;
2004:     mumps->id.rhs_sparse  = NULL;
2005:   }
2006:   mumps->id.ICNTL(20)   = 1; /* rhs is sparse */
2007:   mumps->id.ICNTL(21)   = 0; /* solution is in assembled centralized format */

2009:   /* solve phase */
2010:   /*-------------*/
2011:   mumps->id.job = JOB_SOLVE;
2012:   PetscMUMPS_c(&mumps->id);
2013:   if (mumps->id.INFOG(1) < 0)
2014:     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));

2016:   if (!mumps->myid) {
2017:     MatSeqAIJRestoreArray(Bt,&aa);
2018:     MatRestoreRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&done);
2019:   }
2020:   return(0);
2021: }

2023: /*@
2024:   MatMumpsGetInverse - Get user-specified set of entries in inverse of A

2026:    Logically Collective on Mat

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

2032:   Output Parameter:
2033: . spRHS - requested entries of inverse of A

2035:    Level: beginner

2037:    References:
2038: .      MUMPS Users' Guide

2040: .seealso: MatGetFactor(), MatCreateTranspose()
2041: @*/
2042: PetscErrorCode MatMumpsGetInverse(Mat F,Mat spRHS)
2043: {

2048:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2049:   PetscUseMethod(F,"MatMumpsGetInverse_C",(Mat,Mat),(F,spRHS));
2050:   return(0);
2051: }

2053: /*@
2054:   MatMumpsGetInfo - Get MUMPS parameter INFO()

2056:    Logically Collective on Mat

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

2062:   Output Parameter:
2063: .  ival - value of MUMPS INFO(icntl)

2065:    Level: beginner

2067:    References:
2068: .      MUMPS Users' Guide

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

2078:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2080:   PetscUseMethod(F,"MatMumpsGetInfo_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));
2081:   return(0);
2082: }

2084: /*@
2085:   MatMumpsGetInfog - Get MUMPS parameter INFOG()

2087:    Logically Collective on Mat

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

2093:   Output Parameter:
2094: .  ival - value of MUMPS INFOG(icntl)

2096:    Level: beginner

2098:    References:
2099: .      MUMPS Users' Guide

2101: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2102: @*/
2103: PetscErrorCode MatMumpsGetInfog(Mat F,PetscInt icntl,PetscInt *ival)
2104: {

2109:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2111:   PetscUseMethod(F,"MatMumpsGetInfog_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));
2112:   return(0);
2113: }

2115: /*@
2116:   MatMumpsGetRinfo - Get MUMPS parameter RINFO()

2118:    Logically Collective on Mat

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

2124:   Output Parameter:
2125: .  val - value of MUMPS RINFO(icntl)

2127:    Level: beginner

2129:    References:
2130: .       MUMPS Users' Guide

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

2140:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2142:   PetscUseMethod(F,"MatMumpsGetRinfo_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));
2143:   return(0);
2144: }

2146: /*@
2147:   MatMumpsGetRinfog - Get MUMPS parameter RINFOG()

2149:    Logically Collective on Mat

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

2155:   Output Parameter:
2156: .  val - value of MUMPS RINFOG(icntl)

2158:    Level: beginner

2160:    References:
2161: .      MUMPS Users' Guide

2163: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2164: @*/
2165: PetscErrorCode MatMumpsGetRinfog(Mat F,PetscInt icntl,PetscReal *val)
2166: {

2171:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2173:   PetscUseMethod(F,"MatMumpsGetRinfog_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));
2174:   return(0);
2175: }

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

2181:   Works with MATAIJ and MATSBAIJ matrices

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

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

2187:   Options Database Keys:
2188: +  -mat_mumps_icntl_1 - ICNTL(1): output stream for error messages 
2189: .  -mat_mumps_icntl_2 - ICNTL(2): output stream for diagnostic printing, statistics, and warning 
2190: .  -mat_mumps_icntl_3 -  ICNTL(3): output stream for global information, collected on the host 
2191: .  -mat_mumps_icntl_4 -  ICNTL(4): level of printing (0 to 4) 
2192: .  -mat_mumps_icntl_6 - ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7) 
2193: .  -mat_mumps_icntl_7 - ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis 
2194: .  -mat_mumps_icntl_8  - ICNTL(8): scaling strategy (-2 to 8 or 77) 
2195: .  -mat_mumps_icntl_10  - ICNTL(10): max num of refinements 
2196: .  -mat_mumps_icntl_11  - ICNTL(11): statistics related to an error analysis (via -ksp_view) 
2197: .  -mat_mumps_icntl_12  - ICNTL(12): an ordering strategy for symmetric matrices (0 to 3) 
2198: .  -mat_mumps_icntl_13  - ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting 
2199: .  -mat_mumps_icntl_14  - ICNTL(14): percentage increase in the estimated working space 
2200: .  -mat_mumps_icntl_19  - ICNTL(19): computes the Schur complement 
2201: .  -mat_mumps_icntl_22  - ICNTL(22): in-core/out-of-core factorization and solve (0 or 1) 
2202: .  -mat_mumps_icntl_23  - ICNTL(23): max size of the working memory (MB) that can allocate per processor 
2203: .  -mat_mumps_icntl_24  - ICNTL(24): detection of null pivot rows (0 or 1) 
2204: .  -mat_mumps_icntl_25  - ICNTL(25): compute a solution of a deficient matrix and a null space basis 
2205: .  -mat_mumps_icntl_26  - ICNTL(26): drives the solution phase if a Schur complement matrix 
2206: .  -mat_mumps_icntl_28  - ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering 
2207: .  -mat_mumps_icntl_29 - ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis 
2208: .  -mat_mumps_icntl_30 - ICNTL(30): compute user-specified set of entries in inv(A) 
2209: .  -mat_mumps_icntl_31 - ICNTL(31): indicates which factors may be discarded during factorization 
2210: .  -mat_mumps_icntl_33 - ICNTL(33): compute determinant 
2211: .  -mat_mumps_cntl_1  - CNTL(1): relative pivoting threshold 
2212: .  -mat_mumps_cntl_2  -  CNTL(2): stopping criterion of refinement 
2213: .  -mat_mumps_cntl_3 - CNTL(3): absolute pivoting threshold 
2214: .  -mat_mumps_cntl_4 - CNTL(4): value for static pivoting 
2215: -  -mat_mumps_cntl_5 - CNTL(5): fixation for null pivots 

2217:   Level: beginner

2219:     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 
2220: $          KSPGetPC(ksp,&pc);
2221: $          PCFactorGetMatrix(pc,&mat);
2222: $          MatMumpsGetInfo(mat,....);
2223: $          MatMumpsGetInfog(mat,....); etc.
2224:            Or you can run with -ksp_error_if_not_converged and the program will be stopped and the information printed in the error message.

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

2228: M*/

2230: static PetscErrorCode MatFactorGetSolverType_mumps(Mat A,MatSolverType *type)
2231: {
2233:   *type = MATSOLVERMUMPS;
2234:   return(0);
2235: }

2237: /* MatGetFactor for Seq and MPI AIJ matrices */
2238: static PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
2239: {
2240:   Mat            B;
2242:   Mat_MUMPS      *mumps;
2243:   PetscBool      isSeqAIJ;

2246:   /* Create the factorization matrix */
2247:   PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);
2248:   MatCreate(PetscObjectComm((PetscObject)A),&B);
2249:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
2250:   PetscStrallocpy("mumps",&((PetscObject)B)->type_name);
2251:   MatSetUp(B);

2253:   PetscNewLog(B,&mumps);

2255:   B->ops->view        = MatView_MUMPS;
2256:   B->ops->getinfo     = MatGetInfo_MUMPS;

2258:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);
2259:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);
2260:   PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);
2261:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);
2262:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);
2263:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);
2264:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);
2265:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);
2266:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);
2267:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);
2268:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);
2269:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverse_C",MatMumpsGetInverse_MUMPS);

2271:   if (ftype == MAT_FACTOR_LU) {
2272:     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
2273:     B->factortype            = MAT_FACTOR_LU;
2274:     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
2275:     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
2276:     mumps->sym = 0;
2277:   } else {
2278:     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2279:     B->factortype                  = MAT_FACTOR_CHOLESKY;
2280:     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
2281:     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
2282: #if defined(PETSC_USE_COMPLEX)
2283:     mumps->sym = 2;
2284: #else
2285:     if (A->spd_set && A->spd) mumps->sym = 1;
2286:     else                      mumps->sym = 2;
2287: #endif
2288:   }

2290:   /* set solvertype */
2291:   PetscFree(B->solvertype);
2292:   PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);

2294:   B->ops->destroy = MatDestroy_MUMPS;
2295:   B->data         = (void*)mumps;

2297:   PetscInitializeMUMPS(A,mumps);

2299:   *F = B;
2300:   return(0);
2301: }

2303: /* MatGetFactor for Seq and MPI SBAIJ matrices */
2304: static PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
2305: {
2306:   Mat            B;
2308:   Mat_MUMPS      *mumps;
2309:   PetscBool      isSeqSBAIJ;

2312:   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
2313:   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");
2314:   PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);
2315:   /* Create the factorization matrix */
2316:   MatCreate(PetscObjectComm((PetscObject)A),&B);
2317:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
2318:   PetscStrallocpy("mumps",&((PetscObject)B)->type_name);
2319:   MatSetUp(B);

2321:   PetscNewLog(B,&mumps);
2322:   if (isSeqSBAIJ) {
2323:     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
2324:   } else {
2325:     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
2326:   }

2328:   B->ops->getinfo                = MatGetInfo_External;
2329:   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2330:   B->ops->view                   = MatView_MUMPS;

2332:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);
2333:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);
2334:   PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);
2335:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);
2336:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);
2337:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);
2338:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);
2339:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);
2340:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);
2341:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);
2342:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);
2343:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverse_C",MatMumpsGetInverse_MUMPS);

2345:   B->factortype = MAT_FACTOR_CHOLESKY;
2346: #if defined(PETSC_USE_COMPLEX)
2347:   mumps->sym = 2;
2348: #else
2349:   if (A->spd_set && A->spd) mumps->sym = 1;
2350:   else                      mumps->sym = 2;
2351: #endif

2353:   /* set solvertype */
2354:   PetscFree(B->solvertype);
2355:   PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);

2357:   B->ops->destroy = MatDestroy_MUMPS;
2358:   B->data         = (void*)mumps;

2360:   PetscInitializeMUMPS(A,mumps);

2362:   *F = B;
2363:   return(0);
2364: }

2366: static PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
2367: {
2368:   Mat            B;
2370:   Mat_MUMPS      *mumps;
2371:   PetscBool      isSeqBAIJ;

2374:   /* Create the factorization matrix */
2375:   PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);
2376:   MatCreate(PetscObjectComm((PetscObject)A),&B);
2377:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
2378:   PetscStrallocpy("mumps",&((PetscObject)B)->type_name);
2379:   MatSetUp(B);

2381:   PetscNewLog(B,&mumps);
2382:   if (ftype == MAT_FACTOR_LU) {
2383:     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
2384:     B->factortype            = MAT_FACTOR_LU;
2385:     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
2386:     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
2387:     mumps->sym = 0;
2388:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n");

2390:   B->ops->getinfo     = MatGetInfo_External;
2391:   B->ops->view        = MatView_MUMPS;

2393:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);
2394:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);
2395:   PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);
2396:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);
2397:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);
2398:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);
2399:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);
2400:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);
2401:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);
2402:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);
2403:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);
2404:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverse_C",MatMumpsGetInverse_MUMPS);

2406:   /* set solvertype */
2407:   PetscFree(B->solvertype);
2408:   PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);

2410:   B->ops->destroy = MatDestroy_MUMPS;
2411:   B->data         = (void*)mumps;

2413:   PetscInitializeMUMPS(A,mumps);

2415:   *F = B;
2416:   return(0);
2417: }

2419: /* MatGetFactor for Seq and MPI SELL matrices */
2420: static PetscErrorCode MatGetFactor_sell_mumps(Mat A,MatFactorType ftype,Mat *F)
2421: {
2422:   Mat            B;
2424:   Mat_MUMPS      *mumps;
2425:   PetscBool      isSeqSELL;

2428:   /* Create the factorization matrix */
2429:   PetscObjectTypeCompare((PetscObject)A,MATSEQSELL,&isSeqSELL);
2430:   MatCreate(PetscObjectComm((PetscObject)A),&B);
2431:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
2432:   PetscStrallocpy("mumps",&((PetscObject)B)->type_name);
2433:   MatSetUp(B);

2435:   PetscNewLog(B,&mumps);

2437:   B->ops->view        = MatView_MUMPS;
2438:   B->ops->getinfo     = MatGetInfo_MUMPS;

2440:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);
2441:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);
2442:   PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);
2443:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);
2444:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);
2445:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);
2446:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);
2447:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);
2448:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);
2449:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);
2450:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);

2452:   if (ftype == MAT_FACTOR_LU) {
2453:     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
2454:     B->factortype            = MAT_FACTOR_LU;
2455:     if (isSeqSELL) mumps->ConvertToTriples = MatConvertToTriples_seqsell_seqaij;
2456:     else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"To be implemented");
2457:     mumps->sym = 0;
2458:   } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"To be implemented");

2460:   /* set solvertype */
2461:   PetscFree(B->solvertype);
2462:   PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);

2464:   B->ops->destroy = MatDestroy_MUMPS;
2465:   B->data         = (void*)mumps;

2467:   PetscInitializeMUMPS(A,mumps);

2469:   *F = B;
2470:   return(0);
2471: }

2473: PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_MUMPS(void)
2474: {

2478:   MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);
2479:   MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);
2480:   MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);
2481:   MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);
2482:   MatSolverTypeRegister(MATSOLVERMUMPS,MATMPISBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);
2483:   MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);
2484:   MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);
2485:   MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);
2486:   MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);
2487:   MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);
2488:   MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQSELL,MAT_FACTOR_LU,MatGetFactor_sell_mumps);
2489:   return(0);
2490: }