Actual source code: mumps.c

petsc-master 2019-08-18
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 MUMPS_c cmumps_c
 35: #else
 36: #define MUMPS_c zmumps_c
 37: #endif
 38: #else
 39: #if defined(PETSC_USE_REAL_SINGLE)
 40: #define MUMPS_c smumps_c
 41: #else
 42: #define MUMPS_c dmumps_c
 43: #endif
 44: #endif

 46: /* if using PETSc OpenMP support, we only call MUMPS on master ranks. Before/after the call, we change/restore CPUs the master ranks can run on */
 47: #if defined(PETSC_HAVE_OPENMP_SUPPORT)
 48: #define PetscMUMPS_c(mumps) \
 49:   do { \
 50:     if (mumps->use_petsc_omp_support) { \
 51:       if (mumps->is_omp_master) { \
 52:         PetscOmpCtrlOmpRegionOnMasterBegin(mumps->omp_ctrl); \
 53:         MUMPS_c(&mumps->id); \
 54:         PetscOmpCtrlOmpRegionOnMasterEnd(mumps->omp_ctrl); \
 55:       } \
 56:       PetscOmpCtrlBarrier(mumps->omp_ctrl); \
 57:     } else { \
 58:       MUMPS_c(&mumps->id); \
 59:     } \
 60:   } while(0)
 61: #else
 62: #define PetscMUMPS_c(mumps) \
 63:   do { MUMPS_c(&mumps->id); } while (0)
 64: #endif

 66: /* declare MumpsScalar */
 67: #if defined(PETSC_USE_COMPLEX)
 68: #if defined(PETSC_USE_REAL_SINGLE)
 69: #define MumpsScalar mumps_complex
 70: #else
 71: #define MumpsScalar mumps_double_complex
 72: #endif
 73: #else
 74: #define MumpsScalar PetscScalar
 75: #endif

 77: /* macros s.t. indices match MUMPS documentation */
 78: #define ICNTL(I) icntl[(I)-1]
 79: #define CNTL(I) cntl[(I)-1]
 80: #define INFOG(I) infog[(I)-1]
 81: #define INFO(I) info[(I)-1]
 82: #define RINFOG(I) rinfog[(I)-1]
 83: #define RINFO(I) rinfo[(I)-1]

 85: typedef struct {
 86: #if defined(PETSC_USE_COMPLEX)
 87: #if defined(PETSC_USE_REAL_SINGLE)
 88:   CMUMPS_STRUC_C id;
 89: #else
 90:   ZMUMPS_STRUC_C id;
 91: #endif
 92: #else
 93: #if defined(PETSC_USE_REAL_SINGLE)
 94:   SMUMPS_STRUC_C id;
 95: #else
 96:   DMUMPS_STRUC_C id;
 97: #endif
 98: #endif

100:   MatStructure matstruc;
101:   PetscMPIInt  myid,petsc_size;
102:   PetscInt     *irn,*jcn,nz,sym;
103:   PetscScalar  *val;
104:   MPI_Comm     mumps_comm;
105:   PetscInt     ICNTL9_pre;           /* check if ICNTL(9) is changed from previous MatSolve */
106:   VecScatter   scat_rhs, scat_sol;   /* used by MatSolve() */
107:   Vec          b_seq,x_seq;
108:   PetscInt     ninfo,*info;          /* display INFO */
109:   PetscInt     sizeredrhs;
110:   PetscScalar  *schur_sol;
111:   PetscInt     schur_sizesol;

113:   PetscBool    use_petsc_omp_support;
114:   PetscOmpCtrl omp_ctrl;             /* an OpenMP controler that blocked processes will release their CPU (MPI_Barrier does not have this guarantee) */
115:   MPI_Comm     petsc_comm,omp_comm;  /* petsc_comm is petsc matrix's comm */
116:   PetscMPIInt  mpinz;                /* on master rank, nz = sum(mpinz) over omp_comm; on other ranks, mpinz = nz*/
117:   PetscMPIInt  omp_comm_size;
118:   PetscBool    is_omp_master;        /* is this rank the master of omp_comm */
119:   PetscMPIInt  *recvcount,*displs;

121:   PetscErrorCode (*ConvertToTriples)(Mat, int, MatReuse, int*, int**, int**, PetscScalar**);
122: } Mat_MUMPS;

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

126: static PetscErrorCode MatMumpsResetSchur_Private(Mat_MUMPS* mumps)
127: {

131:   PetscFree2(mumps->id.listvar_schur,mumps->id.schur);
132:   PetscFree(mumps->id.redrhs);
133:   PetscFree(mumps->schur_sol);
134:   mumps->id.size_schur = 0;
135:   mumps->id.schur_lld  = 0;
136:   mumps->id.ICNTL(19)  = 0;
137:   return(0);
138: }

140: /* solve with rhs in mumps->id.redrhs and return in the same location */
141: static PetscErrorCode MatMumpsSolveSchur_Private(Mat F)
142: {
143:   Mat_MUMPS            *mumps=(Mat_MUMPS*)F->data;
144:   Mat                  S,B,X;
145:   MatFactorSchurStatus schurstatus;
146:   PetscInt             sizesol;
147:   PetscErrorCode       ierr;

150:   MatFactorFactorizeSchurComplement(F);
151:   MatFactorGetSchurComplement(F,&S,&schurstatus);
152:   MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.nrhs,(PetscScalar*)mumps->id.redrhs,&B);
153:   switch (schurstatus) {
154:   case MAT_FACTOR_SCHUR_FACTORED:
155:     MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.nrhs,(PetscScalar*)mumps->id.redrhs,&X);
156:     if (!mumps->id.ICNTL(9)) { /* transpose solve */
157:       MatMatSolveTranspose(S,B,X);
158:     } else {
159:       MatMatSolve(S,B,X);
160:     }
161:     break;
162:   case MAT_FACTOR_SCHUR_INVERTED:
163:     sizesol = mumps->id.nrhs*mumps->id.size_schur;
164:     if (!mumps->schur_sol || sizesol > mumps->schur_sizesol) {
165:       PetscFree(mumps->schur_sol);
166:       PetscMalloc1(sizesol,&mumps->schur_sol);
167:       mumps->schur_sizesol = sizesol;
168:     }
169:     MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.nrhs,mumps->schur_sol,&X);
170:     if (!mumps->id.ICNTL(9)) { /* transpose solve */
171:       MatTransposeMatMult(S,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);
172:     } else {
173:       MatMatMult(S,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);
174:     }
175:     MatCopy(X,B,SAME_NONZERO_PATTERN);
176:     break;
177:   default:
178:     SETERRQ1(PetscObjectComm((PetscObject)F),PETSC_ERR_SUP,"Unhandled MatFactorSchurStatus %D",F->schur_status);
179:     break;
180:   }
181:   MatFactorRestoreSchurComplement(F,&S,schurstatus);
182:   MatDestroy(&B);
183:   MatDestroy(&X);
184:   return(0);
185: }

187: static PetscErrorCode MatMumpsHandleSchur_Private(Mat F, PetscBool expansion)
188: {
189:   Mat_MUMPS     *mumps=(Mat_MUMPS*)F->data;

193:   if (!mumps->id.ICNTL(19)) { /* do nothing when Schur complement has not been computed */
194:     return(0);
195:   }
196:   if (!expansion) { /* prepare for the condensation step */
197:     PetscInt sizeredrhs = mumps->id.nrhs*mumps->id.size_schur;
198:     /* allocate MUMPS internal array to store reduced right-hand sides */
199:     if (!mumps->id.redrhs || sizeredrhs > mumps->sizeredrhs) {
200:       PetscFree(mumps->id.redrhs);
201:       mumps->id.lredrhs = mumps->id.size_schur;
202:       PetscMalloc1(mumps->id.nrhs*mumps->id.lredrhs,&mumps->id.redrhs);
203:       mumps->sizeredrhs = mumps->id.nrhs*mumps->id.lredrhs;
204:     }
205:     mumps->id.ICNTL(26) = 1; /* condensation phase */
206:   } else { /* prepare for the expansion step */
207:     /* solve Schur complement (this has to be done by the MUMPS user, so basically us) */
208:     MatMumpsSolveSchur_Private(F);
209:     mumps->id.ICNTL(26) = 2; /* expansion phase */
210:     PetscMUMPS_c(mumps);
211:     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));
212:     /* restore defaults */
213:     mumps->id.ICNTL(26) = -1;
214:     /* free MUMPS internal array for redrhs if we have solved for multiple rhs in order to save memory space */
215:     if (mumps->id.nrhs > 1) {
216:       PetscFree(mumps->id.redrhs);
217:       mumps->id.lredrhs = 0;
218:       mumps->sizeredrhs = 0;
219:     }
220:   }
221:   return(0);
222: }

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

227:   input:
228:     A       - matrix in aij,baij or sbaij format
229:     shift   - 0: C style output triple; 1: Fortran style output triple.
230:     reuse   - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple
231:               MAT_REUSE_MATRIX:   only the values in v array are updated
232:   output:
233:     nnz     - dim of r, c, and v (number of local nonzero entries of A)
234:     r, c, v - row and col index, matrix values (matrix triples)

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

240:  */

242: PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
243: {
244:   const PetscInt *ai,*aj,*ajj,M=A->rmap->n;
245:   PetscInt       nz,rnz,i,j;
247:   PetscInt       *row,*col;
248:   Mat_SeqAIJ     *aa=(Mat_SeqAIJ*)A->data;

251:   *v=aa->a;
252:   if (reuse == MAT_INITIAL_MATRIX) {
253:     nz   = aa->nz;
254:     ai   = aa->i;
255:     aj   = aa->j;
256:     *nnz = nz;
257:     PetscMalloc1(2*nz, &row);
258:     col  = row + nz;

260:     nz = 0;
261:     for (i=0; i<M; i++) {
262:       rnz = ai[i+1] - ai[i];
263:       ajj = aj + ai[i];
264:       for (j=0; j<rnz; j++) {
265:         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
266:       }
267:     }
268:     *r = row; *c = col;
269:   }
270:   return(0);
271: }

273: PetscErrorCode MatConvertToTriples_seqsell_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
274: {
275:   Mat_SeqSELL *a=(Mat_SeqSELL*)A->data;
276:   PetscInt    *ptr;

279:   *v = a->val;
280:   if (reuse == MAT_INITIAL_MATRIX) {
281:     PetscInt       nz,i,j,row;

284:     nz   = a->sliidx[a->totalslices];
285:     *nnz = nz;
286:     PetscMalloc1(2*nz, &ptr);
287:     *r   = ptr;
288:     *c   = ptr + nz;

290:     for (i=0; i<a->totalslices; i++) {
291:       for (j=a->sliidx[i],row=0; j<a->sliidx[i+1]; j++,row=((row+1)&0x07)) {
292:         *ptr++ = 8*i + row + shift;
293:       }
294:     }
295:     for (i=0;i<nz;i++) *ptr++ = a->colidx[i] + shift;
296:   }
297:   return(0);
298: }

300: PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
301: {
302:   Mat_SeqBAIJ    *aa=(Mat_SeqBAIJ*)A->data;
303:   const PetscInt *ai,*aj,*ajj,bs2 = aa->bs2;
304:   PetscInt       bs,M,nz,idx=0,rnz,i,j,k,m;
306:   PetscInt       *row,*col;

309:   MatGetBlockSize(A,&bs);
310:   M = A->rmap->N/bs;
311:   *v = aa->a;
312:   if (reuse == MAT_INITIAL_MATRIX) {
313:     ai   = aa->i; aj = aa->j;
314:     nz   = bs2*aa->nz;
315:     *nnz = nz;
316:     PetscMalloc1(2*nz, &row);
317:     col  = row + nz;

319:     for (i=0; i<M; i++) {
320:       ajj = aj + ai[i];
321:       rnz = ai[i+1] - ai[i];
322:       for (k=0; k<rnz; k++) {
323:         for (j=0; j<bs; j++) {
324:           for (m=0; m<bs; m++) {
325:             row[idx]   = i*bs + m + shift;
326:             col[idx++] = bs*(ajj[k]) + j + shift;
327:           }
328:         }
329:       }
330:     }
331:     *r = row; *c = col;
332:   }
333:   return(0);
334: }

336: PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r,int **c,PetscScalar **v)
337: {
338:   const PetscInt *ai, *aj,*ajj;
339:   PetscInt        nz,rnz,i,j,k,m,bs;
340:   PetscErrorCode  ierr;
341:   PetscInt        *row,*col;
342:   PetscScalar     *val;
343:   Mat_SeqSBAIJ    *aa=(Mat_SeqSBAIJ*)A->data;
344:   const PetscInt  bs2=aa->bs2,mbs=aa->mbs;
345: #if defined(PETSC_USE_COMPLEX)
346:   PetscBool      hermitian;
347: #endif

350: #if defined(PETSC_USE_COMPLEX)
351:   MatGetOption(A,MAT_HERMITIAN,&hermitian);
352:   if (hermitian) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MUMPS does not support Hermitian symmetric matrices for Choleksy");
353: #endif
354:   ai   = aa->i;
355:   aj   = aa->j;
356:   MatGetBlockSize(A,&bs);
357:   if (reuse == MAT_INITIAL_MATRIX) {
358:     nz   = aa->nz;
359:     PetscMalloc((2*bs2*nz*sizeof(PetscInt)+(bs>1?bs2*nz*sizeof(PetscScalar):0)), &row);
360:     col  = row + bs2*nz;
361:     if (bs>1)
362:       val = (PetscScalar*)(col + bs2*nz);
363:     else
364:       val = aa->a;

366:     *r = row; *c = col; *v = val;
367:   } else {
368:     row = *r; col = *c; val = *v;
369:   }

371:   nz = 0;
372:   if (bs>1) {
373:     for (i=0; i<mbs; i++) {
374:       rnz = ai[i+1] - ai[i];
375:       ajj = aj + ai[i];
376:       for (j=0; j<rnz; j++) {
377:         for (k=0; k<bs; k++) {
378:           for (m=0; m<bs; m++) {
379:             if (ajj[j]>i || k>=m) {
380:               if (reuse == MAT_INITIAL_MATRIX) {
381:                 row[nz] = i*bs + m + shift;
382:                 col[nz] = ajj[j]*bs + k + shift;
383:               }
384:               val[nz++] = aa->a[(ai[i]+j)*bs2 + m + k*bs];
385:             }
386:           }
387:         }
388:       }
389:     }
390:   } else if (reuse == MAT_INITIAL_MATRIX) {
391:     for (i=0; i<mbs; i++) {
392:       rnz = ai[i+1] - ai[i];
393:       ajj = aj + ai[i];
394:       for (j=0; j<rnz; j++) {
395:         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
396:       }
397:     }
398:     if (nz != aa->nz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Different numbers of nonzeros %D != %D",nz,aa->nz);
399:   }
400:   *nnz = nz;
401:   return(0);
402: }

404: PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
405: {
406:   const PetscInt    *ai,*aj,*ajj,*adiag,M=A->rmap->n;
407:   PetscInt          nz,rnz,i,j;
408:   const PetscScalar *av,*v1;
409:   PetscScalar       *val;
410:   PetscErrorCode    ierr;
411:   PetscInt          *row,*col;
412:   Mat_SeqAIJ        *aa=(Mat_SeqAIJ*)A->data;
413:   PetscBool         missing;
414: #if defined(PETSC_USE_COMPLEX)
415:   PetscBool         hermitian;
416: #endif

419: #if defined(PETSC_USE_COMPLEX)
420:   MatGetOption(A,MAT_HERMITIAN,&hermitian);
421:   if (hermitian) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MUMPS does not support Hermitian symmetric matrices for Choleksy");
422: #endif
423:   ai    = aa->i; aj = aa->j; av = aa->a;
424:   adiag = aa->diag;
425:   MatMissingDiagonal_SeqAIJ(A,&missing,&i);
426:   if (reuse == MAT_INITIAL_MATRIX) {
427:     /* count nz in the upper triangular part of A */
428:     nz = 0;
429:     if (missing) {
430:       for (i=0; i<M; i++) {
431:         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
432:           for (j=ai[i];j<ai[i+1];j++) {
433:             if (aj[j] < i) continue;
434:             nz++;
435:           }
436:         } else {
437:           nz += ai[i+1] - adiag[i];
438:         }
439:       }
440:     } else {
441:       for (i=0; i<M; i++) nz += ai[i+1] - adiag[i];
442:     }
443:     *nnz = nz;

445:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
446:     col  = row + nz;
447:     val  = (PetscScalar*)(col + nz);

449:     nz = 0;
450:     if (missing) {
451:       for (i=0; i<M; i++) {
452:         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
453:           for (j=ai[i];j<ai[i+1];j++) {
454:             if (aj[j] < i) continue;
455:             row[nz] = i+shift;
456:             col[nz] = aj[j]+shift;
457:             val[nz] = av[j];
458:             nz++;
459:           }
460:         } else {
461:           rnz = ai[i+1] - adiag[i];
462:           ajj = aj + adiag[i];
463:           v1  = av + adiag[i];
464:           for (j=0; j<rnz; j++) {
465:             row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
466:           }
467:         }
468:       }
469:     } else {
470:       for (i=0; i<M; i++) {
471:         rnz = ai[i+1] - adiag[i];
472:         ajj = aj + adiag[i];
473:         v1  = av + adiag[i];
474:         for (j=0; j<rnz; j++) {
475:           row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
476:         }
477:       }
478:     }
479:     *r = row; *c = col; *v = val;
480:   } else {
481:     nz = 0; val = *v;
482:     if (missing) {
483:       for (i=0; i <M; i++) {
484:         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
485:           for (j=ai[i];j<ai[i+1];j++) {
486:             if (aj[j] < i) continue;
487:             val[nz++] = av[j];
488:           }
489:         } else {
490:           rnz = ai[i+1] - adiag[i];
491:           v1  = av + adiag[i];
492:           for (j=0; j<rnz; j++) {
493:             val[nz++] = v1[j];
494:           }
495:         }
496:       }
497:     } else {
498:       for (i=0; i <M; i++) {
499:         rnz = ai[i+1] - adiag[i];
500:         v1  = av + adiag[i];
501:         for (j=0; j<rnz; j++) {
502:           val[nz++] = v1[j];
503:         }
504:       }
505:     }
506:   }
507:   return(0);
508: }

510: PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r,int **c, PetscScalar **v)
511: {
512:   const PetscInt    *ai,*aj,*bi,*bj,*garray,*ajj,*bjj;
513:   PetscErrorCode    ierr;
514:   PetscInt          rstart,nz,bs,i,j,k,m,jj,irow,countA,countB;
515:   PetscInt          *row,*col;
516:   const PetscScalar *av,*bv,*v1,*v2;
517:   PetscScalar       *val;
518:   Mat_MPISBAIJ      *mat = (Mat_MPISBAIJ*)A->data;
519:   Mat_SeqSBAIJ      *aa  = (Mat_SeqSBAIJ*)(mat->A)->data;
520:   Mat_SeqBAIJ       *bb  = (Mat_SeqBAIJ*)(mat->B)->data;
521:   const PetscInt    bs2=aa->bs2,mbs=aa->mbs;
522: #if defined(PETSC_USE_COMPLEX)
523:   PetscBool         hermitian;
524: #endif

527: #if defined(PETSC_USE_COMPLEX)
528:   MatGetOption(A,MAT_HERMITIAN,&hermitian);
529:   if (hermitian) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MUMPS does not support Hermitian symmetric matrices for Choleksy");
530: #endif
531:   MatGetBlockSize(A,&bs);
532:   rstart = A->rmap->rstart;
533:   ai = aa->i;
534:   aj = aa->j;
535:   bi = bb->i;
536:   bj = bb->j;
537:   av = aa->a;
538:   bv = bb->a;

540:   garray = mat->garray;

542:   if (reuse == MAT_INITIAL_MATRIX) {
543:     nz   = aa->nz + bb->nz;
544:     PetscMalloc((2*bs2*nz*sizeof(PetscInt)+bs2*nz*sizeof(PetscScalar)), &row);
545:     col  = row + bs2*nz;
546:     val  = (PetscScalar*)(col + bs2*nz);

548:     *r = row; *c = col; *v = val;
549:   } else {
550:     row = *r; col = *c; val = *v;
551:   }

553:   jj = 0; irow = rstart;
554:   for (i=0; i<mbs; i++) {
555:     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
556:     countA = ai[i+1] - ai[i];
557:     countB = bi[i+1] - bi[i];
558:     bjj    = bj + bi[i];
559:     v1     = av + ai[i]*bs2;
560:     v2     = bv + bi[i]*bs2;

562:     if (bs>1) {
563:       /* A-part */
564:       for (j=0; j<countA; j++) {
565:         for (k=0; k<bs; k++) {
566:           for (m=0; m<bs; m++) {
567:             if (rstart + ajj[j]*bs>irow || k>=m) {
568:               if (reuse == MAT_INITIAL_MATRIX) {
569:                 row[jj] = irow + m + shift; col[jj] = rstart + ajj[j]*bs + k + shift;
570:               }
571:               val[jj++] = v1[j*bs2 + m + k*bs];
572:             }
573:           }
574:         }
575:       }

577:       /* B-part */
578:       for (j=0; j < countB; j++) {
579:         for (k=0; k<bs; k++) {
580:           for (m=0; m<bs; m++) {
581:             if (reuse == MAT_INITIAL_MATRIX) {
582:               row[jj] = irow + m + shift; col[jj] = garray[bjj[j]]*bs + k + shift;
583:             }
584:             val[jj++] = v2[j*bs2 + m + k*bs];
585:           }
586:         }
587:       }
588:     } else {
589:       /* A-part */
590:       for (j=0; j<countA; j++) {
591:         if (reuse == MAT_INITIAL_MATRIX) {
592:           row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
593:         }
594:         val[jj++] = v1[j];
595:       }

597:       /* B-part */
598:       for (j=0; j < countB; j++) {
599:         if (reuse == MAT_INITIAL_MATRIX) {
600:           row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
601:         }
602:         val[jj++] = v2[j];
603:       }
604:     }
605:     irow+=bs;
606:   }
607:   *nnz = jj;
608:   return(0);
609: }

611: PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
612: {
613:   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
614:   PetscErrorCode    ierr;
615:   PetscInt          rstart,nz,i,j,jj,irow,countA,countB;
616:   PetscInt          *row,*col;
617:   const PetscScalar *av, *bv,*v1,*v2;
618:   PetscScalar       *val;
619:   Mat_MPIAIJ        *mat = (Mat_MPIAIJ*)A->data;
620:   Mat_SeqAIJ        *aa  = (Mat_SeqAIJ*)(mat->A)->data;
621:   Mat_SeqAIJ        *bb  = (Mat_SeqAIJ*)(mat->B)->data;

624:   rstart = A->rmap->rstart;
625:   ai = aa->i;
626:   aj = aa->j;
627:   bi = bb->i;
628:   bj = bb->j;
629:   av = aa->a;
630:   bv = bb->a;

632:   garray = mat->garray;

634:   if (reuse == MAT_INITIAL_MATRIX) {
635:     nz   = aa->nz + bb->nz;
636:     *nnz = nz;
637:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
638:     col  = row + nz;
639:     val  = (PetscScalar*)(col + nz);

641:     *r = row; *c = col; *v = val;
642:   } else {
643:     row = *r; col = *c; val = *v;
644:   }

646:   jj = 0; irow = rstart;
647:   for (i=0; i<m; i++) {
648:     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
649:     countA = ai[i+1] - ai[i];
650:     countB = bi[i+1] - bi[i];
651:     bjj    = bj + bi[i];
652:     v1     = av + ai[i];
653:     v2     = bv + bi[i];

655:     /* A-part */
656:     for (j=0; j<countA; j++) {
657:       if (reuse == MAT_INITIAL_MATRIX) {
658:         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
659:       }
660:       val[jj++] = v1[j];
661:     }

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

675: PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
676: {
677:   Mat_MPIBAIJ       *mat    = (Mat_MPIBAIJ*)A->data;
678:   Mat_SeqBAIJ       *aa     = (Mat_SeqBAIJ*)(mat->A)->data;
679:   Mat_SeqBAIJ       *bb     = (Mat_SeqBAIJ*)(mat->B)->data;
680:   const PetscInt    *ai     = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj;
681:   const PetscInt    *garray = mat->garray,mbs=mat->mbs,rstart=A->rmap->rstart;
682:   const PetscInt    bs2=mat->bs2;
683:   PetscErrorCode    ierr;
684:   PetscInt          bs,nz,i,j,k,n,jj,irow,countA,countB,idx;
685:   PetscInt          *row,*col;
686:   const PetscScalar *av=aa->a, *bv=bb->a,*v1,*v2;
687:   PetscScalar       *val;

690:   MatGetBlockSize(A,&bs);
691:   if (reuse == MAT_INITIAL_MATRIX) {
692:     nz   = bs2*(aa->nz + bb->nz);
693:     *nnz = nz;
694:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
695:     col  = row + nz;
696:     val  = (PetscScalar*)(col + nz);

698:     *r = row; *c = col; *v = val;
699:   } else {
700:     row = *r; col = *c; val = *v;
701:   }

703:   jj = 0; irow = rstart;
704:   for (i=0; i<mbs; i++) {
705:     countA = ai[i+1] - ai[i];
706:     countB = bi[i+1] - bi[i];
707:     ajj    = aj + ai[i];
708:     bjj    = bj + bi[i];
709:     v1     = av + bs2*ai[i];
710:     v2     = bv + bs2*bi[i];

712:     idx = 0;
713:     /* A-part */
714:     for (k=0; k<countA; k++) {
715:       for (j=0; j<bs; j++) {
716:         for (n=0; n<bs; n++) {
717:           if (reuse == MAT_INITIAL_MATRIX) {
718:             row[jj] = irow + n + shift;
719:             col[jj] = rstart + bs*ajj[k] + j + shift;
720:           }
721:           val[jj++] = v1[idx++];
722:         }
723:       }
724:     }

726:     idx = 0;
727:     /* B-part */
728:     for (k=0; k<countB; k++) {
729:       for (j=0; j<bs; j++) {
730:         for (n=0; n<bs; n++) {
731:           if (reuse == MAT_INITIAL_MATRIX) {
732:             row[jj] = irow + n + shift;
733:             col[jj] = bs*garray[bjj[k]] + j + shift;
734:           }
735:           val[jj++] = v2[idx++];
736:         }
737:       }
738:     }
739:     irow += bs;
740:   }
741:   return(0);
742: }

744: PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
745: {
746:   const PetscInt    *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
747:   PetscErrorCode    ierr;
748:   PetscInt          rstart,nz,nza,nzb,i,j,jj,irow,countA,countB;
749:   PetscInt          *row,*col;
750:   const PetscScalar *av, *bv,*v1,*v2;
751:   PetscScalar       *val;
752:   Mat_MPIAIJ        *mat =  (Mat_MPIAIJ*)A->data;
753:   Mat_SeqAIJ        *aa  =(Mat_SeqAIJ*)(mat->A)->data;
754:   Mat_SeqAIJ        *bb  =(Mat_SeqAIJ*)(mat->B)->data;
755: #if defined(PETSC_USE_COMPLEX)
756:   PetscBool         hermitian;
757: #endif

760: #if defined(PETSC_USE_COMPLEX)
761:   MatGetOption(A,MAT_HERMITIAN,&hermitian);
762:   if (hermitian) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MUMPS does not support Hermitian symmetric matrices for Choleksy");
763: #endif
764:   ai     = aa->i;
765:   aj     = aa->j;
766:   adiag  = aa->diag;
767:   bi     = bb->i;
768:   bj     = bb->j;
769:   garray = mat->garray;
770:   av     = aa->a;
771:   bv     = bb->a;

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

775:   if (reuse == MAT_INITIAL_MATRIX) {
776:     nza = 0;    /* num of upper triangular entries in mat->A, including diagonals */
777:     nzb = 0;    /* num of upper triangular entries in mat->B */
778:     for (i=0; i<m; i++) {
779:       nza   += (ai[i+1] - adiag[i]);
780:       countB = bi[i+1] - bi[i];
781:       bjj    = bj + bi[i];
782:       for (j=0; j<countB; j++) {
783:         if (garray[bjj[j]] > rstart) nzb++;
784:       }
785:     }

787:     nz   = nza + nzb; /* total nz of upper triangular part of mat */
788:     *nnz = nz;
789:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
790:     col  = row + nz;
791:     val  = (PetscScalar*)(col + nz);

793:     *r = row; *c = col; *v = val;
794:   } else {
795:     row = *r; col = *c; val = *v;
796:   }

798:   jj = 0; irow = rstart;
799:   for (i=0; i<m; i++) {
800:     ajj    = aj + adiag[i];                 /* ptr to the beginning of the diagonal of this row */
801:     v1     = av + adiag[i];
802:     countA = ai[i+1] - adiag[i];
803:     countB = bi[i+1] - bi[i];
804:     bjj    = bj + bi[i];
805:     v2     = bv + bi[i];

807:     /* A-part */
808:     for (j=0; j<countA; j++) {
809:       if (reuse == MAT_INITIAL_MATRIX) {
810:         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
811:       }
812:       val[jj++] = v1[j];
813:     }

815:     /* B-part */
816:     for (j=0; j < countB; j++) {
817:       if (garray[bjj[j]] > rstart) {
818:         if (reuse == MAT_INITIAL_MATRIX) {
819:           row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
820:         }
821:         val[jj++] = v2[j];
822:       }
823:     }
824:     irow++;
825:   }
826:   return(0);
827: }

829: PetscErrorCode MatDestroy_MUMPS(Mat A)
830: {
831:   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->data;

835:   PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);
836:   VecScatterDestroy(&mumps->scat_rhs);
837:   VecScatterDestroy(&mumps->scat_sol);
838:   VecDestroy(&mumps->b_seq);
839:   VecDestroy(&mumps->x_seq);
840:   PetscFree(mumps->id.perm_in);
841:   PetscFree(mumps->irn);
842:   PetscFree(mumps->info);
843:   MatMumpsResetSchur_Private(mumps);
844:   mumps->id.job = JOB_END;
845:   PetscMUMPS_c(mumps);
846: #if defined(PETSC_HAVE_OPENMP_SUPPORT)
847:   if (mumps->use_petsc_omp_support) { PetscOmpCtrlDestroy(&mumps->omp_ctrl); }
848: #endif
849:   PetscFree2(mumps->recvcount,mumps->displs);
850:   PetscFree(A->data);

852:   /* clear composed functions */
853:   PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverType_C",NULL);
854:   PetscObjectComposeFunction((PetscObject)A,"MatFactorSetSchurIS_C",NULL);
855:   PetscObjectComposeFunction((PetscObject)A,"MatFactorCreateSchurComplement_C",NULL);
856:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetIcntl_C",NULL);
857:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetIcntl_C",NULL);
858:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetCntl_C",NULL);
859:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetCntl_C",NULL);
860:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfo_C",NULL);
861:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfog_C",NULL);
862:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfo_C",NULL);
863:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfog_C",NULL);
864:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInverse_C",NULL);
865:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInverseTranspose_C",NULL);
866:   return(0);
867: }

869: PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x)
870: {
871:   Mat_MUMPS        *mumps=(Mat_MUMPS*)A->data;
872:   PetscScalar      *array;
873:   Vec              b_seq;
874:   IS               is_iden,is_petsc;
875:   PetscErrorCode   ierr;
876:   PetscInt         i;
877:   PetscBool        second_solve = PETSC_FALSE;
878:   static PetscBool cite1 = PETSC_FALSE,cite2 = PETSC_FALSE;

881:   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);
882:   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);

884:   if (A->factorerrortype) {
885:     PetscInfo2(A,"MatSolve is called with singular matrix factor, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
886:     VecSetInf(x);
887:     return(0);
888:   }

890:   mumps->id.ICNTL(20)= 0; /* dense RHS */
891:   mumps->id.nrhs     = 1;
892:   b_seq          = mumps->b_seq;
893:   if (mumps->petsc_size > 1) {
894:     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
895:     VecScatterBegin(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);
896:     VecScatterEnd(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);
897:     if (!mumps->myid) {VecGetArray(b_seq,&array);}
898:   } else {  /* petsc_size == 1 */
899:     VecCopy(b,x);
900:     VecGetArray(x,&array);
901:   }
902:   if (!mumps->myid) { /* define rhs on the host */
903:     mumps->id.nrhs = 1;
904:     mumps->id.rhs = (MumpsScalar*)array;
905:   }

907:   /*
908:      handle condensation step of Schur complement (if any)
909:      We set by default ICNTL(26) == -1 when Schur indices have been provided by the user.
910:      According to MUMPS (5.0.0) manual, any value should be harmful during the factorization phase
911:      Unless the user provides a valid value for ICNTL(26), MatSolve and MatMatSolve routines solve the full system.
912:      This requires an extra call to PetscMUMPS_c and the computation of the factors for S
913:   */
914:   if (mumps->id.size_schur > 0 && (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2)) {
915:     if (mumps->petsc_size > 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Parallel Schur complements not yet supported from PETSc\n");
916:     second_solve = PETSC_TRUE;
917:     MatMumpsHandleSchur_Private(A,PETSC_FALSE);
918:   }
919:   /* solve phase */
920:   /*-------------*/
921:   mumps->id.job = JOB_SOLVE;
922:   PetscMUMPS_c(mumps);
923:   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));

925:   /* handle expansion step of Schur complement (if any) */
926:   if (second_solve) {
927:     MatMumpsHandleSchur_Private(A,PETSC_TRUE);
928:   }

930:   if (mumps->petsc_size > 1) { /* convert mumps distributed solution to petsc mpi x */
931:     if (mumps->scat_sol && mumps->ICNTL9_pre != mumps->id.ICNTL(9)) {
932:       /* when id.ICNTL(9) changes, the contents of lsol_loc may change (not its size, lsol_loc), recreates scat_sol */
933:       VecScatterDestroy(&mumps->scat_sol);
934:     }
935:     if (!mumps->scat_sol) { /* create scatter scat_sol */
936:       ISCreateStride(PETSC_COMM_SELF,mumps->id.lsol_loc,0,1,&is_iden); /* from */
937:       for (i=0; i<mumps->id.lsol_loc; i++) {
938:         mumps->id.isol_loc[i] -= 1; /* change Fortran style to C style */
939:       }
940:       ISCreateGeneral(PETSC_COMM_SELF,mumps->id.lsol_loc,mumps->id.isol_loc,PETSC_COPY_VALUES,&is_petsc);  /* to */
941:       VecScatterCreate(mumps->x_seq,is_iden,x,is_petsc,&mumps->scat_sol);
942:       ISDestroy(&is_iden);
943:       ISDestroy(&is_petsc);

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

948:     VecScatterBegin(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);
949:     VecScatterEnd(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);
950:   }

952:   if (mumps->petsc_size > 1) {if (!mumps->myid) {VecRestoreArray(b_seq,&array);}}
953:   else {VecRestoreArray(x,&array);}

955:   PetscLogFlops(2.0*mumps->id.RINFO(3));
956:   return(0);
957: }

959: PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x)
960: {
961:   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->data;

965:   mumps->id.ICNTL(9) = 0;
966:   MatSolve_MUMPS(A,b,x);
967:   mumps->id.ICNTL(9) = 1;
968:   return(0);
969: }

971: PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X)
972: {
973:   PetscErrorCode    ierr;
974:   Mat               Bt = NULL;
975:   PetscBool         flg, flgT;
976:   Mat_MUMPS         *mumps=(Mat_MUMPS*)A->data;
977:   PetscInt          i,nrhs,M;
978:   PetscScalar       *array;
979:   const PetscScalar *rbray;
980:   PetscInt          lsol_loc,nlsol_loc,*isol_loc,*idxx,*isol_loc_save,iidx = 0;
981:   PetscScalar       *bray,*sol_loc,*sol_loc_save;
982:   IS                is_to,is_from;
983:   PetscInt          k,proc,j,m,myrstart;
984:   const PetscInt    *rstart;
985:   Vec               v_mpi,b_seq,msol_loc;
986:   VecScatter        scat_rhs,scat_sol;
987:   PetscScalar       *aa;
988:   PetscInt          spnr,*ia,*ja;
989:   Mat_MPIAIJ        *b = NULL;

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

995:   PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
996:   if (flg) { /* dense B */
997:     if (B->rmap->n != X->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix B and X must have same row distribution");
998:     mumps->id.ICNTL(20)= 0; /* dense RHS */
999:   } else { /* sparse B */
1000:     PetscObjectTypeCompare((PetscObject)B,MATTRANSPOSEMAT,&flgT);
1001:     if (flgT) { /* input B is transpose of actural RHS matrix,
1002:                  because mumps requires sparse compressed COLUMN storage! See MatMatTransposeSolve_MUMPS() */
1003:       MatTransposeGetMat(B,&Bt);
1004:     } else SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B must be MATTRANSPOSEMAT matrix");
1005:     mumps->id.ICNTL(20)= 1; /* sparse RHS */
1006:   }

1008:   MatGetSize(B,&M,&nrhs);
1009:   mumps->id.nrhs = nrhs;
1010:   mumps->id.lrhs = M;
1011:   mumps->id.rhs  = NULL;

1013:   if (mumps->petsc_size == 1) {
1014:     PetscScalar *aa;
1015:     PetscInt    spnr,*ia,*ja;
1016:     PetscBool   second_solve = PETSC_FALSE;

1018:     MatDenseGetArray(X,&array);
1019:     mumps->id.rhs = (MumpsScalar*)array;

1021:     if (!Bt) { /* dense B */
1022:       /* copy B to X */
1023:       MatDenseGetArrayRead(B,&rbray);
1024:       PetscArraycpy(array,rbray,M*nrhs);
1025:       MatDenseRestoreArrayRead(B,&rbray);
1026:     } else { /* sparse B */
1027:       MatSeqAIJGetArray(Bt,&aa);
1028:       MatGetRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&flg);
1029:       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot get IJ structure");
1030:       /* mumps requires ia and ja start at 1! */
1031:       mumps->id.irhs_ptr    = ia;
1032:       mumps->id.irhs_sparse = ja;
1033:       mumps->id.nz_rhs      = ia[spnr] - 1;
1034:       mumps->id.rhs_sparse  = (MumpsScalar*)aa;
1035:     }
1036:     /* handle condensation step of Schur complement (if any) */
1037:     if (mumps->id.size_schur > 0 && (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2)) {
1038:       second_solve = PETSC_TRUE;
1039:       MatMumpsHandleSchur_Private(A,PETSC_FALSE);
1040:     }
1041:     /* solve phase */
1042:     /*-------------*/
1043:     mumps->id.job = JOB_SOLVE;
1044:     PetscMUMPS_c(mumps);
1045:     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));

1047:     /* handle expansion step of Schur complement (if any) */
1048:     if (second_solve) {
1049:       MatMumpsHandleSchur_Private(A,PETSC_TRUE);
1050:     }
1051:     if (Bt) { /* sparse B */
1052:       MatSeqAIJRestoreArray(Bt,&aa);
1053:       MatRestoreRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&flg);
1054:       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot restore IJ structure");
1055:     }
1056:     MatDenseRestoreArray(X,&array);
1057:     return(0);
1058:   }

1060:   /*--------- parallel case: MUMPS requires rhs B to be centralized on the host! --------*/
1061:   if (mumps->petsc_size > 1 && mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Parallel Schur complements not yet supported from PETSc\n");

1063:   /* create msol_loc to hold mumps local solution */
1064:   isol_loc_save = mumps->id.isol_loc; /* save it for MatSolve() */
1065:   sol_loc_save  = (PetscScalar*)mumps->id.sol_loc;

1067:   lsol_loc  = mumps->id.lsol_loc;
1068:   nlsol_loc = nrhs*lsol_loc;     /* length of sol_loc */
1069:   PetscMalloc2(nlsol_loc,&sol_loc,lsol_loc,&isol_loc);
1070:   mumps->id.sol_loc  = (MumpsScalar*)sol_loc;
1071:   mumps->id.isol_loc = isol_loc;

1073:   VecCreateSeqWithArray(PETSC_COMM_SELF,1,nlsol_loc,(PetscScalar*)sol_loc,&msol_loc);

1075:   if (!Bt) { /* dense B */
1076:     /* TODO: Because of non-contiguous indices, the created vecscatter scat_rhs is not done in MPI_Gather, resulting in
1077:        very inefficient communication. An optimization is to use VecScatterCreateToZero to gather B to rank 0. Then on rank
1078:        0, re-arrange B into desired order, which is a local operation.
1079:      */

1081:     /* scatter v_mpi to b_seq because MUMPS only supports centralized rhs */
1082:     /* wrap dense rhs matrix B into a vector v_mpi */
1083:     MatGetLocalSize(B,&m,NULL);
1084:     MatDenseGetArray(B,&bray);
1085:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)B),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi);
1086:     MatDenseRestoreArray(B,&bray);

1088:     /* scatter v_mpi to b_seq in proc[0]. MUMPS requires rhs to be centralized on the host! */
1089:     if (!mumps->myid) {
1090:       PetscInt *idx;
1091:       /* idx: maps from k-th index of v_mpi to (i,j)-th global entry of B */
1092:       PetscMalloc1(nrhs*M,&idx);
1093:       MatGetOwnershipRanges(B,&rstart);
1094:       k = 0;
1095:       for (proc=0; proc<mumps->petsc_size; proc++){
1096:         for (j=0; j<nrhs; j++){
1097:           for (i=rstart[proc]; i<rstart[proc+1]; i++) idx[k++] = j*M + i;
1098:         }
1099:       }

1101:       VecCreateSeq(PETSC_COMM_SELF,nrhs*M,&b_seq);
1102:       ISCreateGeneral(PETSC_COMM_SELF,nrhs*M,idx,PETSC_OWN_POINTER,&is_to);
1103:       ISCreateStride(PETSC_COMM_SELF,nrhs*M,0,1,&is_from);
1104:     } else {
1105:       VecCreateSeq(PETSC_COMM_SELF,0,&b_seq);
1106:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_to);
1107:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_from);
1108:     }
1109:     VecScatterCreate(v_mpi,is_from,b_seq,is_to,&scat_rhs);
1110:     VecScatterBegin(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);
1111:     ISDestroy(&is_to);
1112:     ISDestroy(&is_from);
1113:     VecScatterEnd(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);

1115:     if (!mumps->myid) { /* define rhs on the host */
1116:       VecGetArray(b_seq,&bray);
1117:       mumps->id.rhs = (MumpsScalar*)bray;
1118:       VecRestoreArray(b_seq,&bray);
1119:     }

1121:   } else { /* sparse B */
1122:     b = (Mat_MPIAIJ*)Bt->data;

1124:     /* wrap dense X into a vector v_mpi */
1125:     MatGetLocalSize(X,&m,NULL);
1126:     MatDenseGetArray(X,&bray);
1127:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)X),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi);
1128:     MatDenseRestoreArray(X,&bray);

1130:     if (!mumps->myid) {
1131:       MatSeqAIJGetArray(b->A,&aa);
1132:       MatGetRowIJ(b->A,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&flg);
1133:       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot get IJ structure");
1134:       /* mumps requires ia and ja start at 1! */
1135:       mumps->id.irhs_ptr    = ia;
1136:       mumps->id.irhs_sparse = ja;
1137:       mumps->id.nz_rhs      = ia[spnr] - 1;
1138:       mumps->id.rhs_sparse  = (MumpsScalar*)aa;
1139:     } else {
1140:       mumps->id.irhs_ptr    = NULL;
1141:       mumps->id.irhs_sparse = NULL;
1142:       mumps->id.nz_rhs      = 0;
1143:       mumps->id.rhs_sparse  = NULL;
1144:     }
1145:   }

1147:   /* solve phase */
1148:   /*-------------*/
1149:   mumps->id.job = JOB_SOLVE;
1150:   PetscMUMPS_c(mumps);
1151:   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));

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

1157:   /* create scatter scat_sol */
1158:   MatGetOwnershipRanges(X,&rstart);
1159:   /* iidx: index for scatter mumps solution to petsc X */

1161:   ISCreateStride(PETSC_COMM_SELF,nlsol_loc,0,1,&is_from);
1162:   PetscMalloc1(nlsol_loc,&idxx);
1163:   for (i=0; i<lsol_loc; i++) {
1164:     isol_loc[i] -= 1; /* change Fortran style to C style. isol_loc[i+j*lsol_loc] contains x[isol_loc[i]] in j-th vector */

1166:     for (proc=0; proc<mumps->petsc_size; proc++){
1167:       if (isol_loc[i] >= rstart[proc] && isol_loc[i] < rstart[proc+1]) {
1168:         myrstart = rstart[proc];
1169:         k        = isol_loc[i] - myrstart;        /* local index on 1st column of petsc vector X */
1170:         iidx     = k + myrstart*nrhs;             /* maps mumps isol_loc[i] to petsc index in X */
1171:         m        = rstart[proc+1] - rstart[proc]; /* rows of X for this proc */
1172:         break;
1173:       }
1174:     }

1176:     for (j=0; j<nrhs; j++) idxx[i+j*lsol_loc] = iidx + j*m;
1177:   }
1178:   ISCreateGeneral(PETSC_COMM_SELF,nlsol_loc,idxx,PETSC_COPY_VALUES,&is_to);
1179:   VecScatterCreate(msol_loc,is_from,v_mpi,is_to,&scat_sol);
1180:   VecScatterBegin(scat_sol,msol_loc,v_mpi,INSERT_VALUES,SCATTER_FORWARD);
1181:   ISDestroy(&is_from);
1182:   ISDestroy(&is_to);
1183:   VecScatterEnd(scat_sol,msol_loc,v_mpi,INSERT_VALUES,SCATTER_FORWARD);
1184:   MatDenseRestoreArray(X,&array);

1186:   /* free spaces */
1187:   mumps->id.sol_loc  = (MumpsScalar*)sol_loc_save;
1188:   mumps->id.isol_loc = isol_loc_save;

1190:   PetscFree2(sol_loc,isol_loc);
1191:   PetscFree(idxx);
1192:   VecDestroy(&msol_loc);
1193:   VecDestroy(&v_mpi);
1194:   if (Bt) {
1195:     if (!mumps->myid) {
1196:       b = (Mat_MPIAIJ*)Bt->data;
1197:       MatSeqAIJRestoreArray(b->A,&aa);
1198:       MatRestoreRowIJ(b->A,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&flg);
1199:       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot restore IJ structure");
1200:     }
1201:   } else {
1202:     VecDestroy(&b_seq);
1203:     VecScatterDestroy(&scat_rhs);
1204:   }
1205:   VecScatterDestroy(&scat_sol);
1206:   PetscLogFlops(2.0*nrhs*mumps->id.RINFO(3));
1207:   return(0);
1208: }

1210: PetscErrorCode MatMatTransposeSolve_MUMPS(Mat A,Mat Bt,Mat X)
1211: {
1213:   PetscBool      flg;
1214:   Mat            B;

1217:   PetscObjectTypeCompareAny((PetscObject)Bt,&flg,MATSEQAIJ,MATMPIAIJ,NULL);
1218:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)Bt),PETSC_ERR_ARG_WRONG,"Matrix Bt must be MATAIJ matrix");

1220:   /* Create B=Bt^T that uses Bt's data structure */
1221:   MatCreateTranspose(Bt,&B);

1223:   MatMatSolve_MUMPS(A,B,X);
1224:   MatDestroy(&B);
1225:   return(0);
1226: }

1228: #if !defined(PETSC_USE_COMPLEX)
1229: /*
1230:   input:
1231:    F:        numeric factor
1232:   output:
1233:    nneg:     total number of negative pivots
1234:    nzero:    total number of zero pivots
1235:    npos:     (global dimension of F) - nneg - nzero
1236: */
1237: PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
1238: {
1239:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1241:   PetscMPIInt    size;

1244:   MPI_Comm_size(PetscObjectComm((PetscObject)F),&size);
1245:   /* 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 */
1246:   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));

1248:   if (nneg) *nneg = mumps->id.INFOG(12);
1249:   if (nzero || npos) {
1250:     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");
1251:     if (nzero) *nzero = mumps->id.INFOG(28);
1252:     if (npos) *npos   = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28));
1253:   }
1254:   return(0);
1255: }
1256: #endif

1258: PetscErrorCode MatMumpsGatherNonzerosOnMaster(MatReuse reuse,Mat_MUMPS *mumps)
1259: {
1261:   PetscInt       i,nz=0,*irn,*jcn=0;
1262:   PetscScalar    *val=0;
1263:   PetscMPIInt    mpinz,*recvcount=NULL,*displs=NULL;

1266:   if (mumps->omp_comm_size > 1) {
1267:     if (reuse == MAT_INITIAL_MATRIX) {
1268:       /* master first gathers counts of nonzeros to receive */
1269:       if (mumps->is_omp_master) { PetscMalloc2(mumps->omp_comm_size,&recvcount,mumps->omp_comm_size,&displs); }
1270:       PetscMPIIntCast(mumps->nz,&mpinz);
1271:       MPI_Gather(&mpinz,1,MPI_INT,recvcount,1,MPI_INT,0/*root*/,mumps->omp_comm);

1273:       /* master allocates memory to receive nonzeros */
1274:       if (mumps->is_omp_master) {
1275:         displs[0] = 0;
1276:         for (i=1; i<mumps->omp_comm_size; i++) displs[i] = displs[i-1] + recvcount[i-1];
1277:         nz   = displs[mumps->omp_comm_size-1] + recvcount[mumps->omp_comm_size-1];
1278:         PetscMalloc(2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar),&irn);
1279:         jcn  = irn + nz;
1280:         val  = (PetscScalar*)(jcn + nz);
1281:       }

1283:       /* save the gatherv plan */
1284:       mumps->mpinz     = mpinz; /* used as send count */
1285:       mumps->recvcount = recvcount;
1286:       mumps->displs    = displs;

1288:       /* master gathers nonzeros */
1289:       MPI_Gatherv(mumps->irn,mpinz,MPIU_INT,irn,mumps->recvcount,mumps->displs,MPIU_INT,0/*root*/,mumps->omp_comm);
1290:       MPI_Gatherv(mumps->jcn,mpinz,MPIU_INT,jcn,mumps->recvcount,mumps->displs,MPIU_INT,0/*root*/,mumps->omp_comm);
1291:       MPI_Gatherv(mumps->val,mpinz,MPIU_SCALAR,val,mumps->recvcount,mumps->displs,MPIU_SCALAR,0/*root*/,mumps->omp_comm);

1293:       /* master frees its row/col/val and replaces them with bigger arrays */
1294:       if (mumps->is_omp_master) {
1295:         PetscFree(mumps->irn); /* irn/jcn/val are allocated together so free only irn */
1296:         mumps->nz  = nz; /* it is a sum of mpinz over omp_comm */
1297:         mumps->irn = irn;
1298:         mumps->jcn = jcn;
1299:         mumps->val = val;
1300:       }
1301:     } else {
1302:       MPI_Gatherv((mumps->is_omp_master?MPI_IN_PLACE:mumps->val),mumps->mpinz,MPIU_SCALAR,mumps->val,mumps->recvcount,mumps->displs,MPIU_SCALAR,0/*root*/,mumps->omp_comm);
1303:     }
1304:   }
1305:   return(0);
1306: }

1308: PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
1309: {
1310:   Mat_MUMPS      *mumps =(Mat_MUMPS*)(F)->data;
1312:   PetscBool      isMPIAIJ;

1315:   if (mumps->id.INFOG(1) < 0) {
1316:     if (mumps->id.INFOG(1) == -6) {
1317:       PetscInfo2(A,"MatFactorNumeric is called with singular matrix structure, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1318:     }
1319:     PetscInfo2(A,"MatFactorNumeric is called after analysis phase fails, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1320:     return(0);
1321:   }

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

1326:   /* numerical factorization phase */
1327:   /*-------------------------------*/
1328:   mumps->id.job = JOB_FACTNUMERIC;
1329:   if (!mumps->id.ICNTL(18)) { /* A is centralized */
1330:     if (!mumps->myid) {
1331:       mumps->id.a = (MumpsScalar*)mumps->val;
1332:     }
1333:   } else {
1334:     mumps->id.a_loc = (MumpsScalar*)mumps->val;
1335:   }
1336:   PetscMUMPS_c(mumps);
1337:   if (mumps->id.INFOG(1) < 0) {
1338:     if (A->erroriffailure) {
1339:       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));
1340:     } else {
1341:       if (mumps->id.INFOG(1) == -10) { /* numerically singular matrix */
1342:         PetscInfo2(F,"matrix is numerically singular, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1343:         F->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1344:       } else if (mumps->id.INFOG(1) == -13) {
1345:         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));
1346:         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1347:       } else if (mumps->id.INFOG(1) == -8 || mumps->id.INFOG(1) == -9 || (-16 < mumps->id.INFOG(1) && mumps->id.INFOG(1) < -10) ) {
1348:         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));
1349:         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1350:       } else {
1351:         PetscInfo2(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1352:         F->factorerrortype = MAT_FACTOR_OTHER;
1353:       }
1354:     }
1355:   }
1356:   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));

1358:   F->assembled    = PETSC_TRUE;
1359:   mumps->matstruc = SAME_NONZERO_PATTERN;
1360:   if (F->schur) { /* reset Schur status to unfactored */
1361: #if defined(PETSC_HAVE_CUDA)
1362:     F->schur->valid_GPU_matrix = PETSC_OFFLOAD_CPU;
1363: #endif
1364:     if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
1365:       mumps->id.ICNTL(19) = 2;
1366:       MatTranspose(F->schur,MAT_INPLACE_MATRIX,&F->schur);
1367:     }
1368:     MatFactorRestoreSchurComplement(F,NULL,MAT_FACTOR_SCHUR_UNFACTORED);
1369:   }

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

1374:   if (!mumps->is_omp_master) mumps->id.INFO(23) = 0;
1375:   if (mumps->petsc_size > 1) {
1376:     PetscInt    lsol_loc;
1377:     PetscScalar *sol_loc;

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

1381:     /* distributed solution; Create x_seq=sol_loc for repeated use */
1382:     if (mumps->x_seq) {
1383:       VecScatterDestroy(&mumps->scat_sol);
1384:       PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);
1385:       VecDestroy(&mumps->x_seq);
1386:     }
1387:     lsol_loc = mumps->id.INFO(23); /* length of sol_loc */
1388:     PetscMalloc2(lsol_loc,&sol_loc,lsol_loc,&mumps->id.isol_loc);
1389:     mumps->id.lsol_loc = lsol_loc;
1390:     mumps->id.sol_loc = (MumpsScalar*)sol_loc;
1391:     VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq);
1392:   }
1393:   PetscLogFlops(mumps->id.RINFO(2));
1394:   return(0);
1395: }

1397: /* Sets MUMPS options from the options database */
1398: PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A)
1399: {
1400:   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1402:   PetscInt       icntl,info[80],i,ninfo=80;
1403:   PetscBool      flg;

1406:   PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MUMPS Options","Mat");
1407:   PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);
1408:   if (flg) mumps->id.ICNTL(1) = icntl;
1409:   PetscOptionsInt("-mat_mumps_icntl_2","ICNTL(2): output stream for diagnostic printing, statistics, and warning","None",mumps->id.ICNTL(2),&icntl,&flg);
1410:   if (flg) mumps->id.ICNTL(2) = icntl;
1411:   PetscOptionsInt("-mat_mumps_icntl_3","ICNTL(3): output stream for global information, collected on the host","None",mumps->id.ICNTL(3),&icntl,&flg);
1412:   if (flg) mumps->id.ICNTL(3) = icntl;

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

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

1421:   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);
1422:   if (flg) {
1423:     if (icntl== 1 && mumps->petsc_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");
1424:     else mumps->id.ICNTL(7) = icntl;
1425:   }

1427:   PetscOptionsInt("-mat_mumps_icntl_8","ICNTL(8): scaling strategy (-2 to 8 or 77)","None",mumps->id.ICNTL(8),&mumps->id.ICNTL(8),NULL);
1428:   /* 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() */
1429:   PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",mumps->id.ICNTL(10),&mumps->id.ICNTL(10),NULL);
1430:   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);
1431:   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);
1432:   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);
1433:   PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage increase in the estimated working space","None",mumps->id.ICNTL(14),&mumps->id.ICNTL(14),NULL);
1434:   PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): computes the Schur complement","None",mumps->id.ICNTL(19),&mumps->id.ICNTL(19),NULL);
1435:   if (mumps->id.ICNTL(19) <= 0 || mumps->id.ICNTL(19) > 3) { /* reset any schur data (if any) */
1436:     MatDestroy(&F->schur);
1437:     MatMumpsResetSchur_Private(mumps);
1438:   }
1439:   /* 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 */
1440:   /* 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 */

1442:   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);
1443:   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);
1444:   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);
1445:   if (mumps->id.ICNTL(24)) {
1446:     mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */
1447:   }

1449:   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);
1450:   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);
1451:   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);
1452:   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);
1453:   PetscOptionsInt("-mat_mumps_icntl_29","ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis","None",mumps->id.ICNTL(29),&mumps->id.ICNTL(29),NULL);
1454:   /* 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 */
1455:   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);
1456:   /* 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 */
1457:   PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),NULL);
1458:   PetscOptionsInt("-mat_mumps_icntl_35","ICNTL(35): activates Block Low Rank (BLR) based factorization","None",mumps->id.ICNTL(35),&mumps->id.ICNTL(35),NULL);
1459:   PetscOptionsInt("-mat_mumps_icntl_36","ICNTL(36): choice of BLR factorization variant","None",mumps->id.ICNTL(36),&mumps->id.ICNTL(36),NULL);
1460:   PetscOptionsInt("-mat_mumps_icntl_38","ICNTL(38): estimated compression rate of LU factors with BLR","None",mumps->id.ICNTL(38),&mumps->id.ICNTL(38),NULL);

1462:   PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),NULL);
1463:   PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",mumps->id.CNTL(2),&mumps->id.CNTL(2),NULL);
1464:   PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),NULL);
1465:   PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",mumps->id.CNTL(4),&mumps->id.CNTL(4),NULL);
1466:   PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",mumps->id.CNTL(5),&mumps->id.CNTL(5),NULL);
1467:   PetscOptionsReal("-mat_mumps_cntl_7","CNTL(7): dropping parameter used during BLR","None",mumps->id.CNTL(7),&mumps->id.CNTL(7),NULL);

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

1471:   PetscOptionsIntArray("-mat_mumps_view_info","request INFO local to each processor","",info,&ninfo,NULL);
1472:   if (ninfo) {
1473:     if (ninfo > 80) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"number of INFO %d must <= 80\n",ninfo);
1474:     PetscMalloc1(ninfo,&mumps->info);
1475:     mumps->ninfo = ninfo;
1476:     for (i=0; i<ninfo; i++) {
1477:       if (info[i] < 0 || info[i]>80) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index of INFO %d must between 1 and 80\n",ninfo);
1478:       else  mumps->info[i] = info[i];
1479:     }
1480:   }

1482:   PetscOptionsEnd();
1483:   return(0);
1484: }

1486: PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps)
1487: {
1489:   PetscInt       nthreads=0;

1492:   mumps->petsc_comm = PetscObjectComm((PetscObject)A);
1493:   MPI_Comm_size(mumps->petsc_comm,&mumps->petsc_size);
1494:   MPI_Comm_rank(mumps->petsc_comm,&mumps->myid); /* so that code like "if (!myid)" still works even if mumps_comm is different */

1496:   PetscOptionsHasName(NULL,NULL,"-mat_mumps_use_omp_threads",&mumps->use_petsc_omp_support);
1497:   if (mumps->use_petsc_omp_support) nthreads = -1; /* -1 will let PetscOmpCtrlCreate() guess a proper value when user did not supply one */
1498:   PetscOptionsGetInt(NULL,NULL,"-mat_mumps_use_omp_threads",&nthreads,NULL);
1499:   if (mumps->use_petsc_omp_support) {
1500: #if defined(PETSC_HAVE_OPENMP_SUPPORT)
1501:     PetscOmpCtrlCreate(mumps->petsc_comm,nthreads,&mumps->omp_ctrl);
1502:     PetscOmpCtrlGetOmpComms(mumps->omp_ctrl,&mumps->omp_comm,&mumps->mumps_comm,&mumps->is_omp_master);
1503: #else
1504:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP_SYS,"the system does not have PETSc OpenMP support but you added the -mat_mumps_use_omp_threads option. Configure PETSc with --with-openmp --download-hwloc (or --with-hwloc) to enable it, see more in MATSOLVERMUMPS manual\n");
1505: #endif
1506:   } else {
1507:     mumps->omp_comm      = PETSC_COMM_SELF;
1508:     mumps->mumps_comm    = mumps->petsc_comm;
1509:     mumps->is_omp_master = PETSC_TRUE;
1510:   }
1511:   MPI_Comm_size(mumps->omp_comm,&mumps->omp_comm_size);

1513:   mumps->id.comm_fortran = MPI_Comm_c2f(mumps->mumps_comm);
1514:   mumps->id.job = JOB_INIT;
1515:   mumps->id.par = 1;  /* host participates factorizaton and solve */
1516:   mumps->id.sym = mumps->sym;

1518:   PetscMUMPS_c(mumps);

1520:   /* copy MUMPS default control values from master to slaves. Although slaves do not call MUMPS, they may access these values in code.
1521:      For example, ICNTL(9) is initialized to 1 by MUMPS and slaves check ICNTL(9) in MatSolve_MUMPS.
1522:    */
1523:   MPI_Bcast(mumps->id.icntl,60,MPIU_INT, 0,mumps->omp_comm); /* see MUMPS-5.1.2 Manual Section 9 */
1524:   MPI_Bcast(mumps->id.cntl, 15,MPIU_REAL,0,mumps->omp_comm);

1526:   mumps->scat_rhs     = NULL;
1527:   mumps->scat_sol     = NULL;

1529:   /* set PETSc-MUMPS default options - override MUMPS default */
1530:   mumps->id.ICNTL(3) = 0;
1531:   mumps->id.ICNTL(4) = 0;
1532:   if (mumps->petsc_size == 1) {
1533:     mumps->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
1534:   } else {
1535:     mumps->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
1536:     mumps->id.ICNTL(20) = 0;   /* rhs is in dense format */
1537:     mumps->id.ICNTL(21) = 1;   /* distributed solution */
1538:   }

1540:   /* schur */
1541:   mumps->id.size_schur      = 0;
1542:   mumps->id.listvar_schur   = NULL;
1543:   mumps->id.schur           = NULL;
1544:   mumps->sizeredrhs         = 0;
1545:   mumps->schur_sol          = NULL;
1546:   mumps->schur_sizesol      = 0;
1547:   return(0);
1548: }

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

1555:   MPI_Bcast(mumps->id.infog, 80,MPIU_INT, 0,mumps->omp_comm); /* see MUMPS-5.1.2 manual p82 */
1556:   MPI_Bcast(mumps->id.rinfog,20,MPIU_REAL,0,mumps->omp_comm);
1557:   if (mumps->id.INFOG(1) < 0) {
1558:     if (A->erroriffailure) {
1559:       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
1560:     } else {
1561:       if (mumps->id.INFOG(1) == -6) {
1562:         PetscInfo2(F,"matrix is singular in structure, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1563:         F->factorerrortype = MAT_FACTOR_STRUCT_ZEROPIVOT;
1564:       } else if (mumps->id.INFOG(1) == -5 || mumps->id.INFOG(1) == -7) {
1565:         PetscInfo2(F,"problem of workspace, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1566:         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1567:       } else {
1568:         PetscInfo2(F,"Error reported by MUMPS in analysis phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1569:         F->factorerrortype = MAT_FACTOR_OTHER;
1570:       }
1571:     }
1572:   }
1573:   return(0);
1574: }

1576: /* Note Petsc r(=c) permutation is used when mumps->id.ICNTL(7)==1 with centralized assembled matrix input; otherwise r and c are ignored */
1577: PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1578: {
1579:   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1581:   Vec            b;
1582:   const PetscInt M = A->rmap->N;

1585:   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;

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

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

1593:   /* analysis phase */
1594:   /*----------------*/
1595:   mumps->id.job = JOB_FACTSYMBOLIC;
1596:   mumps->id.n   = M;
1597:   switch (mumps->id.ICNTL(18)) {
1598:   case 0:  /* centralized assembled matrix input */
1599:     if (!mumps->myid) {
1600:       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1601:       if (mumps->id.ICNTL(6)>1) {
1602:         mumps->id.a = (MumpsScalar*)mumps->val;
1603:       }
1604:       if (mumps->id.ICNTL(7) == 1) { /* use user-provide matrix ordering - assuming r = c ordering */
1605:         /*
1606:         PetscBool      flag;
1607:         ISEqual(r,c,&flag);
1608:         if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"row_perm != col_perm");
1609:         ISView(r,PETSC_VIEWER_STDOUT_SELF);
1610:          */
1611:         if (!mumps->myid) {
1612:           const PetscInt *idx;
1613:           PetscInt       i,*perm_in;

1615:           PetscMalloc1(M,&perm_in);
1616:           ISGetIndices(r,&idx);

1618:           mumps->id.perm_in = perm_in;
1619:           for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */
1620:           ISRestoreIndices(r,&idx);
1621:         }
1622:       }
1623:     }
1624:     break;
1625:   case 3:  /* distributed assembled matrix input (size>1) */
1626:     mumps->id.nz_loc = mumps->nz;
1627:     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1628:     if (mumps->id.ICNTL(6)>1) {
1629:       mumps->id.a_loc = (MumpsScalar*)mumps->val;
1630:     }
1631:     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1632:     MatCreateVecs(A,NULL,&b);
1633:     VecScatterCreateToZero(b,&mumps->scat_rhs,&mumps->b_seq);
1634:     VecDestroy(&b);
1635:     break;
1636:   }
1637:   PetscMUMPS_c(mumps);
1638:   MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);

1640:   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1641:   F->ops->solve           = MatSolve_MUMPS;
1642:   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1643:   F->ops->matsolve        = MatMatSolve_MUMPS;
1644:   F->ops->mattransposesolve = MatMatTransposeSolve_MUMPS;
1645:   return(0);
1646: }

1648: /* Note the Petsc r and c permutations are ignored */
1649: PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1650: {
1651:   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1653:   Vec            b;
1654:   const PetscInt M = A->rmap->N;

1657:   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;

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

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

1665:   /* analysis phase */
1666:   /*----------------*/
1667:   mumps->id.job = JOB_FACTSYMBOLIC;
1668:   mumps->id.n   = M;
1669:   switch (mumps->id.ICNTL(18)) {
1670:   case 0:  /* centralized assembled matrix input */
1671:     if (!mumps->myid) {
1672:       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1673:       if (mumps->id.ICNTL(6)>1) {
1674:         mumps->id.a = (MumpsScalar*)mumps->val;
1675:       }
1676:     }
1677:     break;
1678:   case 3:  /* distributed assembled matrix input (size>1) */
1679:     mumps->id.nz_loc = mumps->nz;
1680:     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1681:     if (mumps->id.ICNTL(6)>1) {
1682:       mumps->id.a_loc = (MumpsScalar*)mumps->val;
1683:     }
1684:     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1685:     MatCreateVecs(A,NULL,&b);
1686:     VecScatterCreateToZero(b,&mumps->scat_rhs,&mumps->b_seq);
1687:     VecDestroy(&b);
1688:     break;
1689:   }
1690:   PetscMUMPS_c(mumps);
1691:   MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);

1693:   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1694:   F->ops->solve           = MatSolve_MUMPS;
1695:   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1696:   return(0);
1697: }

1699: /* Note the Petsc r permutation and factor info are ignored */
1700: PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
1701: {
1702:   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1704:   Vec            b;
1705:   const PetscInt M = A->rmap->N;

1708:   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;

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

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

1716:   /* analysis phase */
1717:   /*----------------*/
1718:   mumps->id.job = JOB_FACTSYMBOLIC;
1719:   mumps->id.n   = M;
1720:   switch (mumps->id.ICNTL(18)) {
1721:   case 0:  /* centralized assembled matrix input */
1722:     if (!mumps->myid) {
1723:       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1724:       if (mumps->id.ICNTL(6)>1) {
1725:         mumps->id.a = (MumpsScalar*)mumps->val;
1726:       }
1727:     }
1728:     break;
1729:   case 3:  /* distributed assembled matrix input (size>1) */
1730:     mumps->id.nz_loc = mumps->nz;
1731:     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1732:     if (mumps->id.ICNTL(6)>1) {
1733:       mumps->id.a_loc = (MumpsScalar*)mumps->val;
1734:     }
1735:     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1736:     MatCreateVecs(A,NULL,&b);
1737:     VecScatterCreateToZero(b,&mumps->scat_rhs,&mumps->b_seq);
1738:     VecDestroy(&b);
1739:     break;
1740:   }
1741:   PetscMUMPS_c(mumps);
1742:   MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);

1744:   F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
1745:   F->ops->solve                 = MatSolve_MUMPS;
1746:   F->ops->solvetranspose        = MatSolve_MUMPS;
1747:   F->ops->matsolve              = MatMatSolve_MUMPS;
1748:   F->ops->mattransposesolve     = MatMatTransposeSolve_MUMPS;
1749: #if defined(PETSC_USE_COMPLEX)
1750:   F->ops->getinertia = NULL;
1751: #else
1752:   F->ops->getinertia = MatGetInertia_SBAIJMUMPS;
1753: #endif
1754:   return(0);
1755: }

1757: PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
1758: {
1759:   PetscErrorCode    ierr;
1760:   PetscBool         iascii;
1761:   PetscViewerFormat format;
1762:   Mat_MUMPS         *mumps=(Mat_MUMPS*)A->data;

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

1768:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1769:   if (iascii) {
1770:     PetscViewerGetFormat(viewer,&format);
1771:     if (format == PETSC_VIEWER_ASCII_INFO) {
1772:       PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");
1773:       PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                   %d \n",mumps->id.sym);
1774:       PetscViewerASCIIPrintf(viewer,"  PAR (host participation):            %d \n",mumps->id.par);
1775:       PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):         %d \n",mumps->id.ICNTL(1));
1776:       PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg): %d \n",mumps->id.ICNTL(2));
1777:       PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):   %d \n",mumps->id.ICNTL(3));
1778:       PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):        %d \n",mumps->id.ICNTL(4));
1779:       PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):         %d \n",mumps->id.ICNTL(5));
1780:       PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):        %d \n",mumps->id.ICNTL(6));
1781:       PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (sequential matrix ordering):%d \n",mumps->id.ICNTL(7));
1782:       PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scaling strategy):        %d \n",mumps->id.ICNTL(8));
1783:       PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements):  %d \n",mumps->id.ICNTL(10));
1784:       PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):          %d \n",mumps->id.ICNTL(11));
1785:       if (mumps->id.ICNTL(11)>0) {
1786:         PetscViewerASCIIPrintf(viewer,"    RINFOG(4) (inf norm of input mat):        %g\n",mumps->id.RINFOG(4));
1787:         PetscViewerASCIIPrintf(viewer,"    RINFOG(5) (inf norm of solution):         %g\n",mumps->id.RINFOG(5));
1788:         PetscViewerASCIIPrintf(viewer,"    RINFOG(6) (inf norm of residual):         %g\n",mumps->id.RINFOG(6));
1789:         PetscViewerASCIIPrintf(viewer,"    RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8));
1790:         PetscViewerASCIIPrintf(viewer,"    RINFOG(9) (error estimate):               %g \n",mumps->id.RINFOG(9));
1791:         PetscViewerASCIIPrintf(viewer,"    RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11));
1792:       }
1793:       PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",mumps->id.ICNTL(12));
1794:       PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",mumps->id.ICNTL(13));
1795:       PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",mumps->id.ICNTL(14));
1796:       /* ICNTL(15-17) not used */
1797:       PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",mumps->id.ICNTL(18));
1798:       PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Schur complement info):                      %d \n",mumps->id.ICNTL(19));
1799:       PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",mumps->id.ICNTL(20));
1800:       PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d \n",mumps->id.ICNTL(21));
1801:       PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",mumps->id.ICNTL(22));
1802:       PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",mumps->id.ICNTL(23));

1804:       PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",mumps->id.ICNTL(24));
1805:       PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",mumps->id.ICNTL(25));
1806:       PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",mumps->id.ICNTL(26));
1807:       PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",mumps->id.ICNTL(27));
1808:       PetscViewerASCIIPrintf(viewer,"  ICNTL(28) (use parallel or sequential ordering):        %d \n",mumps->id.ICNTL(28));
1809:       PetscViewerASCIIPrintf(viewer,"  ICNTL(29) (parallel ordering):                          %d \n",mumps->id.ICNTL(29));

1811:       PetscViewerASCIIPrintf(viewer,"  ICNTL(30) (user-specified set of entries in inv(A)):    %d \n",mumps->id.ICNTL(30));
1812:       PetscViewerASCIIPrintf(viewer,"  ICNTL(31) (factors is discarded in the solve phase):    %d \n",mumps->id.ICNTL(31));
1813:       PetscViewerASCIIPrintf(viewer,"  ICNTL(33) (compute determinant):                        %d \n",mumps->id.ICNTL(33));
1814:       PetscViewerASCIIPrintf(viewer,"  ICNTL(35) (activate BLR based factorization):           %d \n",mumps->id.ICNTL(35));
1815:       PetscViewerASCIIPrintf(viewer,"  ICNTL(36) (choice of BLR factorization variant):        %d \n",mumps->id.ICNTL(36));
1816:       PetscViewerASCIIPrintf(viewer,"  ICNTL(38) (estimated compression rate of LU factors):   %d \n",mumps->id.ICNTL(38));

1818:       PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",mumps->id.CNTL(1));
1819:       PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2));
1820:       PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",mumps->id.CNTL(3));
1821:       PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",mumps->id.CNTL(4));
1822:       PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",mumps->id.CNTL(5));
1823:       PetscViewerASCIIPrintf(viewer,"  CNTL(7) (dropping parameter for BLR):       %g \n",mumps->id.CNTL(7));

1825:       /* infomation local to each processor */
1826:       PetscViewerASCIIPrintf(viewer, "  RINFO(1) (local estimated flops for the elimination after analysis): \n");
1827:       PetscViewerASCIIPushSynchronized(viewer);
1828:       PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %g \n",mumps->myid,mumps->id.RINFO(1));
1829:       PetscViewerFlush(viewer);
1830:       PetscViewerASCIIPrintf(viewer, "  RINFO(2) (local estimated flops for the assembly after factorization): \n");
1831:       PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(2));
1832:       PetscViewerFlush(viewer);
1833:       PetscViewerASCIIPrintf(viewer, "  RINFO(3) (local estimated flops for the elimination after factorization): \n");
1834:       PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(3));
1835:       PetscViewerFlush(viewer);

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

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

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

1849:       if (mumps->ninfo && mumps->ninfo <= 80){
1850:         PetscInt i;
1851:         for (i=0; i<mumps->ninfo; i++){
1852:           PetscViewerASCIIPrintf(viewer, "  INFO(%d): \n",mumps->info[i]);
1853:           PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(mumps->info[i]));
1854:           PetscViewerFlush(viewer);
1855:         }
1856:       }
1857:       PetscViewerASCIIPopSynchronized(viewer);

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

1865:         PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(3));
1866:         PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(4));
1867:         PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",mumps->id.INFOG(5));
1868:         PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",mumps->id.INFOG(6));
1869:         PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively use after analysis): %d \n",mumps->id.INFOG(7));
1870:         PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",mumps->id.INFOG(8));
1871:         PetscViewerASCIIPrintf(viewer,"  INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",mumps->id.INFOG(9));
1872:         PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",mumps->id.INFOG(10));
1873:         PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",mumps->id.INFOG(11));
1874:         PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",mumps->id.INFOG(12));
1875:         PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",mumps->id.INFOG(13));
1876:         PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",mumps->id.INFOG(14));
1877:         PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",mumps->id.INFOG(15));
1878:         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));
1879:         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));
1880:         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));
1881:         PetscViewerASCIIPrintf(viewer,"  INFOG(19) (size of all MUMPS internal data allocated during factorization: sum over all processors): %d \n",mumps->id.INFOG(19));
1882:         PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",mumps->id.INFOG(20));
1883:         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));
1884:         PetscViewerASCIIPrintf(viewer,"  INFOG(22) (size in MB of memory effectively used during factorization - sum over all processors): %d \n",mumps->id.INFOG(22));
1885:         PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",mumps->id.INFOG(23));
1886:         PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",mumps->id.INFOG(24));
1887:         PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",mumps->id.INFOG(25));
1888:         PetscViewerASCIIPrintf(viewer,"  INFOG(28) (after factorization: number of null pivots encountered): %d\n",mumps->id.INFOG(28));
1889:         PetscViewerASCIIPrintf(viewer,"  INFOG(29) (after factorization: effective number of entries in the factors (sum over all processors)): %d\n",mumps->id.INFOG(29));
1890:         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));
1891:         PetscViewerASCIIPrintf(viewer,"  INFOG(32) (after analysis: type of analysis done): %d\n",mumps->id.INFOG(32));
1892:         PetscViewerASCIIPrintf(viewer,"  INFOG(33) (value used for ICNTL(8)): %d\n",mumps->id.INFOG(33));
1893:         PetscViewerASCIIPrintf(viewer,"  INFOG(34) (exponent of the determinant if determinant is requested): %d\n",mumps->id.INFOG(34));
1894:         PetscViewerASCIIPrintf(viewer,"  INFOG(35) (after factorization: number of entries taking into account BLR factor compression - sum over all processors): %d\n",mumps->id.INFOG(35));
1895:         PetscViewerASCIIPrintf(viewer,"  INFOG(36) (after analysis: estimated size of all MUMPS internal data for running BLR in-core - value on the most memory consuming processor): %d \n",mumps->id.INFOG(36));
1896:         PetscViewerASCIIPrintf(viewer,"  INFOG(37) (after analysis: estimated size of all MUMPS internal data for running BLR in-core - sum over all processors): %d \n",mumps->id.INFOG(37));
1897:         PetscViewerASCIIPrintf(viewer,"  INFOG(38) (after analysis: estimated size of all MUMPS internal data for running BLR out-of-core - value on the most memory consuming processor): %d \n",mumps->id.INFOG(38));
1898:         PetscViewerASCIIPrintf(viewer,"  INFOG(39) (after analysis: estimated size of all MUMPS internal data for running BLR out-of-core - sum over all processors): %d \n",mumps->id.INFOG(39));
1899:       }
1900:     }
1901:   }
1902:   return(0);
1903: }

1905: PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
1906: {
1907:   Mat_MUMPS *mumps =(Mat_MUMPS*)A->data;

1910:   info->block_size        = 1.0;
1911:   info->nz_allocated      = mumps->id.INFOG(20);
1912:   info->nz_used           = mumps->id.INFOG(20);
1913:   info->nz_unneeded       = 0.0;
1914:   info->assemblies        = 0.0;
1915:   info->mallocs           = 0.0;
1916:   info->memory            = 0.0;
1917:   info->fill_ratio_given  = 0;
1918:   info->fill_ratio_needed = 0;
1919:   info->factor_mallocs    = 0;
1920:   return(0);
1921: }

1923: /* -------------------------------------------------------------------------------------------*/
1924: PetscErrorCode MatFactorSetSchurIS_MUMPS(Mat F, IS is)
1925: {
1926:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1927:   const PetscInt *idxs;
1928:   PetscInt       size,i;

1932:   ISGetLocalSize(is,&size);
1933:   if (mumps->petsc_size > 1) {
1934:     PetscBool ls,gs; /* gs is false if any rank other than root has non-empty IS */

1936:     ls   = mumps->myid ? (size ? PETSC_FALSE : PETSC_TRUE) : PETSC_TRUE; /* always true on root; false on others if their size != 0 */
1937:     MPI_Allreduce(&ls,&gs,1,MPIU_BOOL,MPI_LAND,mumps->petsc_comm);
1938:     if (!gs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MUMPS distributed parallel Schur complements not yet supported from PETSc\n");
1939:   }
1940:   if (mumps->id.size_schur != size) {
1941:     PetscFree2(mumps->id.listvar_schur,mumps->id.schur);
1942:     mumps->id.size_schur = size;
1943:     mumps->id.schur_lld  = size;
1944:     PetscMalloc2(size,&mumps->id.listvar_schur,size*size,&mumps->id.schur);
1945:   }

1947:   /* Schur complement matrix */
1948:   MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.size_schur,(PetscScalar*)mumps->id.schur,&F->schur);
1949:   if (mumps->sym == 1) {
1950:     MatSetOption(F->schur,MAT_SPD,PETSC_TRUE);
1951:   }

1953:   /* MUMPS expects Fortran style indices */
1954:   ISGetIndices(is,&idxs);
1955:   PetscArraycpy(mumps->id.listvar_schur,idxs,size);
1956:   for (i=0;i<size;i++) mumps->id.listvar_schur[i]++;
1957:   ISRestoreIndices(is,&idxs);
1958:   if (mumps->petsc_size > 1) {
1959:     mumps->id.ICNTL(19) = 1; /* MUMPS returns Schur centralized on the host */
1960:   } else {
1961:     if (F->factortype == MAT_FACTOR_LU) {
1962:       mumps->id.ICNTL(19) = 3; /* MUMPS returns full matrix */
1963:     } else {
1964:       mumps->id.ICNTL(19) = 2; /* MUMPS returns lower triangular part */
1965:     }
1966:   }
1967:   /* set a special value of ICNTL (not handled my MUMPS) to be used in the solve phase by PETSc */
1968:   mumps->id.ICNTL(26) = -1;
1969:   return(0);
1970: }

1972: /* -------------------------------------------------------------------------------------------*/
1973: PetscErrorCode MatFactorCreateSchurComplement_MUMPS(Mat F,Mat* S)
1974: {
1975:   Mat            St;
1976:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1977:   PetscScalar    *array;
1978: #if defined(PETSC_USE_COMPLEX)
1979:   PetscScalar    im = PetscSqrtScalar((PetscScalar)-1.0);
1980: #endif

1984:   if (!mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it");
1985:   MatCreate(PETSC_COMM_SELF,&St);
1986:   MatSetSizes(St,PETSC_DECIDE,PETSC_DECIDE,mumps->id.size_schur,mumps->id.size_schur);
1987:   MatSetType(St,MATDENSE);
1988:   MatSetUp(St);
1989:   MatDenseGetArray(St,&array);
1990:   if (!mumps->sym) { /* MUMPS always return a full matrix */
1991:     if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
1992:       PetscInt i,j,N=mumps->id.size_schur;
1993:       for (i=0;i<N;i++) {
1994:         for (j=0;j<N;j++) {
1995: #if !defined(PETSC_USE_COMPLEX)
1996:           PetscScalar val = mumps->id.schur[i*N+j];
1997: #else
1998:           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
1999: #endif
2000:           array[j*N+i] = val;
2001:         }
2002:       }
2003:     } else { /* stored by columns */
2004:       PetscArraycpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur);
2005:     }
2006:   } else { /* either full or lower-triangular (not packed) */
2007:     if (mumps->id.ICNTL(19) == 2) { /* lower triangular stored by columns */
2008:       PetscInt i,j,N=mumps->id.size_schur;
2009:       for (i=0;i<N;i++) {
2010:         for (j=i;j<N;j++) {
2011: #if !defined(PETSC_USE_COMPLEX)
2012:           PetscScalar val = mumps->id.schur[i*N+j];
2013: #else
2014:           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
2015: #endif
2016:           array[i*N+j] = val;
2017:           array[j*N+i] = val;
2018:         }
2019:       }
2020:     } else if (mumps->id.ICNTL(19) == 3) { /* full matrix */
2021:       PetscArraycpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur);
2022:     } else { /* ICNTL(19) == 1 lower triangular stored by rows */
2023:       PetscInt i,j,N=mumps->id.size_schur;
2024:       for (i=0;i<N;i++) {
2025:         for (j=0;j<i+1;j++) {
2026: #if !defined(PETSC_USE_COMPLEX)
2027:           PetscScalar val = mumps->id.schur[i*N+j];
2028: #else
2029:           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
2030: #endif
2031:           array[i*N+j] = val;
2032:           array[j*N+i] = val;
2033:         }
2034:       }
2035:     }
2036:   }
2037:   MatDenseRestoreArray(St,&array);
2038:   *S   = St;
2039:   return(0);
2040: }

2042: /* -------------------------------------------------------------------------------------------*/
2043: PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival)
2044: {
2045:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

2048:   mumps->id.ICNTL(icntl) = ival;
2049:   return(0);
2050: }

2052: PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt *ival)
2053: {
2054:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

2057:   *ival = mumps->id.ICNTL(icntl);
2058:   return(0);
2059: }

2061: /*@
2062:   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()

2064:    Logically Collective on Mat

2066:    Input Parameters:
2067: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2068: .  icntl - index of MUMPS parameter array ICNTL()
2069: -  ival - value of MUMPS ICNTL(icntl)

2071:   Options Database:
2072: .   -mat_mumps_icntl_<icntl> <ival>

2074:    Level: beginner

2076:    References:
2077: .     MUMPS Users' Guide

2079: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2080:  @*/
2081: PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival)
2082: {

2087:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2090:   PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));
2091:   return(0);
2092: }

2094: /*@
2095:   MatMumpsGetIcntl - Get MUMPS parameter ICNTL()

2097:    Logically Collective on Mat

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

2103:   Output Parameter:
2104: .  ival - value of MUMPS ICNTL(icntl)

2106:    Level: beginner

2108:    References:
2109: .     MUMPS Users' Guide

2111: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2112: @*/
2113: PetscErrorCode MatMumpsGetIcntl(Mat F,PetscInt icntl,PetscInt *ival)
2114: {

2119:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2122:   PetscUseMethod(F,"MatMumpsGetIcntl_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));
2123:   return(0);
2124: }

2126: /* -------------------------------------------------------------------------------------------*/
2127: PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal val)
2128: {
2129:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

2132:   mumps->id.CNTL(icntl) = val;
2133:   return(0);
2134: }

2136: PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal *val)
2137: {
2138:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

2141:   *val = mumps->id.CNTL(icntl);
2142:   return(0);
2143: }

2145: /*@
2146:   MatMumpsSetCntl - Set MUMPS parameter CNTL()

2148:    Logically Collective on Mat

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

2155:   Options Database:
2156: .   -mat_mumps_cntl_<icntl> <val>

2158:    Level: beginner

2160:    References:
2161: .     MUMPS Users' Guide

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

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

2178: /*@
2179:   MatMumpsGetCntl - Get MUMPS parameter CNTL()

2181:    Logically Collective on Mat

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

2187:   Output Parameter:
2188: .  val - value of MUMPS CNTL(icntl)

2190:    Level: beginner

2192:    References:
2193: .      MUMPS Users' Guide

2195: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2196: @*/
2197: PetscErrorCode MatMumpsGetCntl(Mat F,PetscInt icntl,PetscReal *val)
2198: {

2203:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2206:   PetscUseMethod(F,"MatMumpsGetCntl_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));
2207:   return(0);
2208: }

2210: PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F,PetscInt icntl,PetscInt *info)
2211: {
2212:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

2215:   *info = mumps->id.INFO(icntl);
2216:   return(0);
2217: }

2219: PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F,PetscInt icntl,PetscInt *infog)
2220: {
2221:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

2224:   *infog = mumps->id.INFOG(icntl);
2225:   return(0);
2226: }

2228: PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfo)
2229: {
2230:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

2233:   *rinfo = mumps->id.RINFO(icntl);
2234:   return(0);
2235: }

2237: PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfog)
2238: {
2239:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

2242:   *rinfog = mumps->id.RINFOG(icntl);
2243:   return(0);
2244: }

2246: PetscErrorCode MatMumpsGetInverse_MUMPS(Mat F,Mat spRHS)
2247: {
2249:   Mat            Bt = NULL,Btseq = NULL;
2250:   PetscBool      flg;
2251:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
2252:   PetscScalar    *aa;
2253:   PetscInt       spnr,*ia,*ja;

2257:   PetscObjectTypeCompare((PetscObject)spRHS,MATTRANSPOSEMAT,&flg);
2258:   if (flg) {
2259:     MatTransposeGetMat(spRHS,&Bt);
2260:   } else SETERRQ(PetscObjectComm((PetscObject)spRHS),PETSC_ERR_ARG_WRONG,"Matrix spRHS must be type MATTRANSPOSEMAT matrix");

2262:   MatMumpsSetIcntl(F,30,1);

2264:   if (mumps->petsc_size > 1) {
2265:     Mat_MPIAIJ *b = (Mat_MPIAIJ*)Bt->data;
2266:     Btseq = b->A;
2267:   } else {
2268:     Btseq = Bt;
2269:   }

2271:   if (!mumps->myid) {
2272:     MatSeqAIJGetArray(Btseq,&aa);
2273:     MatGetRowIJ(Btseq,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&flg);
2274:     if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot get IJ structure");

2276:     mumps->id.irhs_ptr    = ia;
2277:     mumps->id.irhs_sparse = ja;
2278:     mumps->id.nz_rhs      = ia[spnr] - 1;
2279:     mumps->id.rhs_sparse  = (MumpsScalar*)aa;
2280:   } else {
2281:     mumps->id.irhs_ptr    = NULL;
2282:     mumps->id.irhs_sparse = NULL;
2283:     mumps->id.nz_rhs      = 0;
2284:     mumps->id.rhs_sparse  = NULL;
2285:   }
2286:   mumps->id.ICNTL(20)   = 1; /* rhs is sparse */
2287:   mumps->id.ICNTL(21)   = 0; /* solution is in assembled centralized format */

2289:   /* solve phase */
2290:   /*-------------*/
2291:   mumps->id.job = JOB_SOLVE;
2292:   PetscMUMPS_c(mumps);
2293:   if (mumps->id.INFOG(1) < 0)
2294:     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));

2296:   if (!mumps->myid) {
2297:     MatSeqAIJRestoreArray(Btseq,&aa);
2298:     MatRestoreRowIJ(Btseq,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&flg);
2299:     if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot get IJ structure");
2300:   }
2301:   return(0);
2302: }

2304: /*@
2305:   MatMumpsGetInverse - Get user-specified set of entries in inverse of A

2307:    Logically Collective on Mat

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

2313:   Output Parameter:
2314: . spRHS - requested entries of inverse of A

2316:    Level: beginner

2318:    References:
2319: .      MUMPS Users' Guide

2321: .seealso: MatGetFactor(), MatCreateTranspose()
2322: @*/
2323: PetscErrorCode MatMumpsGetInverse(Mat F,Mat spRHS)
2324: {

2329:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2330:   PetscUseMethod(F,"MatMumpsGetInverse_C",(Mat,Mat),(F,spRHS));
2331:   return(0);
2332: }

2334: PetscErrorCode MatMumpsGetInverseTranspose_MUMPS(Mat F,Mat spRHST)
2335: {
2337:   Mat            spRHS;

2340:   MatCreateTranspose(spRHST,&spRHS);
2341:   MatMumpsGetInverse_MUMPS(F,spRHS);
2342:   MatDestroy(&spRHS);
2343:   return(0);
2344: }

2346: /*@
2347:   MatMumpsGetInverseTranspose - Get user-specified set of entries in inverse of matrix A^T

2349:    Logically Collective on Mat

2351:    Input Parameters:
2352: +  F - the factored matrix of A obtained by calling MatGetFactor() from PETSc-MUMPS interface
2353: -  spRHST - sequential sparse matrix in MATAIJ format holding specified indices of A^T in processor[0]

2355:   Output Parameter:
2356: . spRHST - requested entries of inverse of A^T

2358:    Level: beginner

2360:    References:
2361: .      MUMPS Users' Guide

2363: .seealso: MatGetFactor(), MatCreateTranspose(), MatMumpsGetInverse()
2364: @*/
2365: PetscErrorCode MatMumpsGetInverseTranspose(Mat F,Mat spRHST)
2366: {
2368:   PetscBool      flg;

2372:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2373:   PetscObjectTypeCompareAny((PetscObject)spRHST,&flg,MATSEQAIJ,MATMPIAIJ,NULL);
2374:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)spRHST),PETSC_ERR_ARG_WRONG,"Matrix spRHST must be MATAIJ matrix");

2376:   PetscUseMethod(F,"MatMumpsGetInverseTranspose_C",(Mat,Mat),(F,spRHST));
2377:   return(0);
2378: }

2380: /*@
2381:   MatMumpsGetInfo - Get MUMPS parameter INFO()

2383:    Logically Collective on Mat

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

2389:   Output Parameter:
2390: .  ival - value of MUMPS INFO(icntl)

2392:    Level: beginner

2394:    References:
2395: .      MUMPS Users' Guide

2397: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2398: @*/
2399: PetscErrorCode MatMumpsGetInfo(Mat F,PetscInt icntl,PetscInt *ival)
2400: {

2405:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2407:   PetscUseMethod(F,"MatMumpsGetInfo_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));
2408:   return(0);
2409: }

2411: /*@
2412:   MatMumpsGetInfog - Get MUMPS parameter INFOG()

2414:    Logically Collective on Mat

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

2420:   Output Parameter:
2421: .  ival - value of MUMPS INFOG(icntl)

2423:    Level: beginner

2425:    References:
2426: .      MUMPS Users' Guide

2428: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2429: @*/
2430: PetscErrorCode MatMumpsGetInfog(Mat F,PetscInt icntl,PetscInt *ival)
2431: {

2436:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2438:   PetscUseMethod(F,"MatMumpsGetInfog_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));
2439:   return(0);
2440: }

2442: /*@
2443:   MatMumpsGetRinfo - Get MUMPS parameter RINFO()

2445:    Logically Collective on Mat

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

2451:   Output Parameter:
2452: .  val - value of MUMPS RINFO(icntl)

2454:    Level: beginner

2456:    References:
2457: .       MUMPS Users' Guide

2459: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2460: @*/
2461: PetscErrorCode MatMumpsGetRinfo(Mat F,PetscInt icntl,PetscReal *val)
2462: {

2467:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2469:   PetscUseMethod(F,"MatMumpsGetRinfo_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));
2470:   return(0);
2471: }

2473: /*@
2474:   MatMumpsGetRinfog - Get MUMPS parameter RINFOG()

2476:    Logically Collective on Mat

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

2482:   Output Parameter:
2483: .  val - value of MUMPS RINFOG(icntl)

2485:    Level: beginner

2487:    References:
2488: .      MUMPS Users' Guide

2490: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2491: @*/
2492: PetscErrorCode MatMumpsGetRinfog(Mat F,PetscInt icntl,PetscReal *val)
2493: {

2498:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2500:   PetscUseMethod(F,"MatMumpsGetRinfog_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));
2501:   return(0);
2502: }

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

2508:   Works with MATAIJ and MATSBAIJ matrices

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

2512:   Use ./configure --with-openmp --download-hwloc (or --with-hwloc) to enable running MUMPS in MPI+OpenMP hybrid mode and non-MUMPS in flat-MPI mode. See details below.

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

2516:   Options Database Keys:
2517: +  -mat_mumps_icntl_1 - ICNTL(1): output stream for error messages
2518: .  -mat_mumps_icntl_2 - ICNTL(2): output stream for diagnostic printing, statistics, and warning
2519: .  -mat_mumps_icntl_3 -  ICNTL(3): output stream for global information, collected on the host
2520: .  -mat_mumps_icntl_4 -  ICNTL(4): level of printing (0 to 4)
2521: .  -mat_mumps_icntl_6 - ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7)
2522: .  -mat_mumps_icntl_7 - ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis
2523: .  -mat_mumps_icntl_8  - ICNTL(8): scaling strategy (-2 to 8 or 77)
2524: .  -mat_mumps_icntl_10  - ICNTL(10): max num of refinements
2525: .  -mat_mumps_icntl_11  - ICNTL(11): statistics related to an error analysis (via -ksp_view)
2526: .  -mat_mumps_icntl_12  - ICNTL(12): an ordering strategy for symmetric matrices (0 to 3)
2527: .  -mat_mumps_icntl_13  - ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting
2528: .  -mat_mumps_icntl_14  - ICNTL(14): percentage increase in the estimated working space
2529: .  -mat_mumps_icntl_19  - ICNTL(19): computes the Schur complement
2530: .  -mat_mumps_icntl_22  - ICNTL(22): in-core/out-of-core factorization and solve (0 or 1)
2531: .  -mat_mumps_icntl_23  - ICNTL(23): max size of the working memory (MB) that can allocate per processor
2532: .  -mat_mumps_icntl_24  - ICNTL(24): detection of null pivot rows (0 or 1)
2533: .  -mat_mumps_icntl_25  - ICNTL(25): compute a solution of a deficient matrix and a null space basis
2534: .  -mat_mumps_icntl_26  - ICNTL(26): drives the solution phase if a Schur complement matrix
2535: .  -mat_mumps_icntl_28  - ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering
2536: .  -mat_mumps_icntl_29 - ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis
2537: .  -mat_mumps_icntl_30 - ICNTL(30): compute user-specified set of entries in inv(A)
2538: .  -mat_mumps_icntl_31 - ICNTL(31): indicates which factors may be discarded during factorization
2539: .  -mat_mumps_icntl_33 - ICNTL(33): compute determinant
2540: .  -mat_mumps_icntl_35 - ICNTL(35): level of activation of BLR (Block Low-Rank) feature
2541: .  -mat_mumps_icntl_36 - ICNTL(36): controls the choice of BLR factorization variant
2542: .  -mat_mumps_icntl_38 - ICNTL(38): sets the estimated compression rate of LU factors with BLR
2543: .  -mat_mumps_cntl_1  - CNTL(1): relative pivoting threshold
2544: .  -mat_mumps_cntl_2  -  CNTL(2): stopping criterion of refinement
2545: .  -mat_mumps_cntl_3 - CNTL(3): absolute pivoting threshold
2546: .  -mat_mumps_cntl_4 - CNTL(4): value for static pivoting
2547: .  -mat_mumps_cntl_5 - CNTL(5): fixation for null pivots
2548: .  -mat_mumps_cntl_7 - CNTL(7): precision of the dropping parameter used during BLR factorization
2549: -  -mat_mumps_use_omp_threads [m] - run MUMPS in MPI+OpenMP hybrid mode as if omp_set_num_threads(m) is called before calling MUMPS.
2550:                                    Default might be the number of cores per CPU package (socket) as reported by hwloc and suggested by the MUMPS manual.

2552:   Level: beginner

2554:     Notes:
2555:     MUMPS Cholesky does not handle (complex) Hermitian matrices http://mumps.enseeiht.fr/doc/userguide_5.2.1.pdf so using it will error if the matrix is Hermitian.

2557:     When a MUMPS factorization fails inside a KSP solve, for example with a KSP_DIVERGED_PC_FAILED, one can find the MUMPS information about the failure by calling
2558: $          KSPGetPC(ksp,&pc);
2559: $          PCFactorGetMatrix(pc,&mat);
2560: $          MatMumpsGetInfo(mat,....);
2561: $          MatMumpsGetInfog(mat,....); etc.
2562:            Or you can run with -ksp_error_if_not_converged and the program will be stopped and the information printed in the error message.

2564:    Two modes to run MUMPS/PETSc with OpenMP

2566: $     Set OMP_NUM_THREADS and run with fewer MPI ranks than cores. For example, if you want to have 16 OpenMP
2567: $     threads per rank, then you may use "export OMP_NUM_THREADS=16 && mpirun -n 4 ./test".

2569: $     -mat_mumps_use_omp_threads [m] and run your code with as many MPI ranks as the number of cores. For example,
2570: $     if a compute node has 32 cores and you run on two nodes, you may use "mpirun -n 64 ./test -mat_mumps_use_omp_threads 16"

2572:    To run MUMPS in MPI+OpenMP hybrid mode (i.e., enable multithreading in MUMPS), but still run the non-MUMPS part
2573:    (i.e., PETSc part) of your code in the so-called flat-MPI (aka pure-MPI) mode, you need to configure PETSc with --with-openmp --download-hwloc
2574:    (or --with-hwloc), and have an MPI that supports MPI-3.0's process shared memory (which is usually available). Since MUMPS calls BLAS
2575:    libraries, to really get performance, you should have multithreaded BLAS libraries such as Intel MKL, AMD ACML, Cray libSci or OpenBLAS
2576:    (PETSc will automatically try to utilized a threaded BLAS if --with-openmp is provided).

2578:    If you run your code through a job submission system, there are caveats in MPI rank mapping. We use MPI_Comm_split_type() to obtain MPI
2579:    processes on each compute node. Listing the processes in rank ascending order, we split processes on a node into consecutive groups of
2580:    size m and create a communicator called omp_comm for each group. Rank 0 in an omp_comm is called the master rank, and others in the omp_comm
2581:    are called slave ranks (or slaves). Only master ranks are seen to MUMPS and slaves are not. We will free CPUs assigned to slaves (might be set
2582:    by CPU binding policies in job scripts) and make the CPUs available to the master so that OMP threads spawned by MUMPS can run on the CPUs.
2583:    In a multi-socket compute node, MPI rank mapping is an issue. Still use the above example and suppose your compute node has two sockets,
2584:    if you interleave MPI ranks on the two sockets, in other words, even ranks are placed on socket 0, and odd ranks are on socket 1, and bind
2585:    MPI ranks to cores, then with -mat_mumps_use_omp_threads 16, a master rank (and threads it spawns) will use half cores in socket 0, and half
2586:    cores in socket 1, that definitely hurts locality. On the other hand, if you map MPI ranks consecutively on the two sockets, then the
2587:    problem will not happen. Therefore, when you use -mat_mumps_use_omp_threads, you need to keep an eye on your MPI rank mapping and CPU binding.
2588:    For example, with the Slurm job scheduler, one can use srun --cpu-bind=verbose -m block:block to map consecutive MPI ranks to sockets and
2589:    examine the mapping result.

2591:    PETSc does not control thread binding in MUMPS. So to get best performance, one still has to set OMP_PROC_BIND and OMP_PLACES in job scripts,
2592:    for example, export OMP_PLACES=threads and export OMP_PROC_BIND=spread. One does not need to export OMP_NUM_THREADS=m in job scripts as PETSc
2593:    calls omp_set_num_threads(m) internally before calling MUMPS.

2595:    References:
2596: +   1. - Heroux, Michael A., R. Brightwell, and Michael M. Wolf. "Bi-modal MPI and MPI+ threads computing on scalable multicore systems." IJHPCA (Submitted) (2011).
2597: -   2. - Gutierrez, Samuel K., et al. "Accommodating Thread-Level Heterogeneity in Coupled Parallel Applications." Parallel and Distributed Processing Symposium (IPDPS), 2017 IEEE International. IEEE, 2017.

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

2601: M*/

2603: static PetscErrorCode MatFactorGetSolverType_mumps(Mat A,MatSolverType *type)
2604: {
2606:   *type = MATSOLVERMUMPS;
2607:   return(0);
2608: }

2610: /* MatGetFactor for Seq and MPI AIJ matrices */
2611: static PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
2612: {
2613:   Mat            B;
2615:   Mat_MUMPS      *mumps;
2616:   PetscBool      isSeqAIJ;

2619:   /* Create the factorization matrix */
2620:   PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);
2621:   MatCreate(PetscObjectComm((PetscObject)A),&B);
2622:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
2623:   PetscStrallocpy("mumps",&((PetscObject)B)->type_name);
2624:   MatSetUp(B);

2626:   PetscNewLog(B,&mumps);

2628:   B->ops->view        = MatView_MUMPS;
2629:   B->ops->getinfo     = MatGetInfo_MUMPS;

2631:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);
2632:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);
2633:   PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);
2634:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);
2635:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);
2636:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);
2637:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);
2638:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);
2639:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);
2640:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);
2641:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);
2642:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverse_C",MatMumpsGetInverse_MUMPS);
2643:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverseTranspose_C",MatMumpsGetInverseTranspose_MUMPS);

2645:   if (ftype == MAT_FACTOR_LU) {
2646:     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
2647:     B->factortype            = MAT_FACTOR_LU;
2648:     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
2649:     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
2650:     mumps->sym = 0;
2651:   } else {
2652:     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2653:     B->factortype                  = MAT_FACTOR_CHOLESKY;
2654:     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
2655:     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
2656: #if defined(PETSC_USE_COMPLEX)
2657:     mumps->sym = 2;
2658: #else
2659:     if (A->spd_set && A->spd) mumps->sym = 1;
2660:     else                      mumps->sym = 2;
2661: #endif
2662:   }

2664:   /* set solvertype */
2665:   PetscFree(B->solvertype);
2666:   PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);

2668:   B->ops->destroy = MatDestroy_MUMPS;
2669:   B->data         = (void*)mumps;

2671:   PetscInitializeMUMPS(A,mumps);

2673:   *F = B;
2674:   return(0);
2675: }

2677: /* MatGetFactor for Seq and MPI SBAIJ matrices */
2678: static PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
2679: {
2680:   Mat            B;
2682:   Mat_MUMPS      *mumps;
2683:   PetscBool      isSeqSBAIJ;

2686:   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
2687:   PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);
2688:   /* Create the factorization matrix */
2689:   MatCreate(PetscObjectComm((PetscObject)A),&B);
2690:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
2691:   PetscStrallocpy("mumps",&((PetscObject)B)->type_name);
2692:   MatSetUp(B);

2694:   PetscNewLog(B,&mumps);
2695:   if (isSeqSBAIJ) {
2696:     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
2697:   } else {
2698:     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
2699:   }

2701:   B->ops->getinfo                = MatGetInfo_External;
2702:   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2703:   B->ops->view                   = MatView_MUMPS;

2705:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);
2706:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);
2707:   PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);
2708:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);
2709:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);
2710:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);
2711:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);
2712:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);
2713:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);
2714:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);
2715:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);
2716:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverse_C",MatMumpsGetInverse_MUMPS);
2717:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverseTranspose_C",MatMumpsGetInverseTranspose_MUMPS);

2719:   B->factortype = MAT_FACTOR_CHOLESKY;
2720: #if defined(PETSC_USE_COMPLEX)
2721:   mumps->sym = 2;
2722: #else
2723:   if (A->spd_set && A->spd) mumps->sym = 1;
2724:   else                      mumps->sym = 2;
2725: #endif

2727:   /* set solvertype */
2728:   PetscFree(B->solvertype);
2729:   PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);

2731:   B->ops->destroy = MatDestroy_MUMPS;
2732:   B->data         = (void*)mumps;

2734:   PetscInitializeMUMPS(A,mumps);

2736:   *F = B;
2737:   return(0);
2738: }

2740: static PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
2741: {
2742:   Mat            B;
2744:   Mat_MUMPS      *mumps;
2745:   PetscBool      isSeqBAIJ;

2748:   /* Create the factorization matrix */
2749:   PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);
2750:   MatCreate(PetscObjectComm((PetscObject)A),&B);
2751:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
2752:   PetscStrallocpy("mumps",&((PetscObject)B)->type_name);
2753:   MatSetUp(B);

2755:   PetscNewLog(B,&mumps);
2756:   if (ftype == MAT_FACTOR_LU) {
2757:     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
2758:     B->factortype            = MAT_FACTOR_LU;
2759:     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
2760:     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
2761:     mumps->sym = 0;
2762:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n");

2764:   B->ops->getinfo     = MatGetInfo_External;
2765:   B->ops->view        = MatView_MUMPS;

2767:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);
2768:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);
2769:   PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);
2770:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);
2771:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);
2772:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);
2773:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);
2774:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);
2775:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);
2776:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);
2777:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);
2778:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverse_C",MatMumpsGetInverse_MUMPS);
2779:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverseTranspose_C",MatMumpsGetInverseTranspose_MUMPS);

2781:   /* set solvertype */
2782:   PetscFree(B->solvertype);
2783:   PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);

2785:   B->ops->destroy = MatDestroy_MUMPS;
2786:   B->data         = (void*)mumps;

2788:   PetscInitializeMUMPS(A,mumps);

2790:   *F = B;
2791:   return(0);
2792: }

2794: /* MatGetFactor for Seq and MPI SELL matrices */
2795: static PetscErrorCode MatGetFactor_sell_mumps(Mat A,MatFactorType ftype,Mat *F)
2796: {
2797:   Mat            B;
2799:   Mat_MUMPS      *mumps;
2800:   PetscBool      isSeqSELL;

2803:   /* Create the factorization matrix */
2804:   PetscObjectTypeCompare((PetscObject)A,MATSEQSELL,&isSeqSELL);
2805:   MatCreate(PetscObjectComm((PetscObject)A),&B);
2806:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
2807:   PetscStrallocpy("mumps",&((PetscObject)B)->type_name);
2808:   MatSetUp(B);

2810:   PetscNewLog(B,&mumps);

2812:   B->ops->view        = MatView_MUMPS;
2813:   B->ops->getinfo     = MatGetInfo_MUMPS;

2815:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);
2816:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);
2817:   PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);
2818:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);
2819:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);
2820:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);
2821:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);
2822:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);
2823:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);
2824:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);
2825:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);

2827:   if (ftype == MAT_FACTOR_LU) {
2828:     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
2829:     B->factortype            = MAT_FACTOR_LU;
2830:     if (isSeqSELL) mumps->ConvertToTriples = MatConvertToTriples_seqsell_seqaij;
2831:     else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"To be implemented");
2832:     mumps->sym = 0;
2833:   } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"To be implemented");

2835:   /* set solvertype */
2836:   PetscFree(B->solvertype);
2837:   PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);

2839:   B->ops->destroy = MatDestroy_MUMPS;
2840:   B->data         = (void*)mumps;

2842:   PetscInitializeMUMPS(A,mumps);

2844:   *F = B;
2845:   return(0);
2846: }

2848: PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_MUMPS(void)
2849: {

2853:   MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);
2854:   MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);
2855:   MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);
2856:   MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);
2857:   MatSolverTypeRegister(MATSOLVERMUMPS,MATMPISBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);
2858:   MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);
2859:   MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);
2860:   MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);
2861:   MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);
2862:   MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);
2863:   MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQSELL,MAT_FACTOR_LU,MatGetFactor_sell_mumps);
2864:   return(0);
2865: }