Actual source code: mumps.c

petsc-master 2020-06-02
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: /* MUMPS uses MUMPS_INT for nonzero indices such as irn/jcn, irn_loc/jcn_loc and uses int64_t for
 47:    number of nonzeros such as nnz, nnz_loc. We typedef MUMPS_INT to PetscMUMPSInt to follow the
 48:    naming convention in PetscMPIInt, PetscBLASInt etc.
 49: */
 50: typedef MUMPS_INT PetscMUMPSInt;

 52: #if defined(INTSIZE64)            /* INTSIZE64 is a macro one used to build MUMPS in full 64-bit mode */
 53: #error "Petsc has not been tested with full 64-bit MUMPS and we choose to error out"
 54: #else
 55: #define MPIU_MUMPSINT             MPI_INT
 56: #define PETSC_MUMPS_INT_MAX       2147483647
 57: #define PETSC_MUMPS_INT_MIN       -2147483648
 58: #endif

 60: /* Cast PetscInt to PetscMUMPSInt. Usually there is no overflow since <a> is row/col indices or some small integers*/
 61: PETSC_STATIC_INLINE PetscErrorCode PetscMUMPSIntCast(PetscInt a,PetscMUMPSInt *b)
 62: {
 64:   if (PetscDefined(USE_64BIT_INDICES) && PetscUnlikelyDebug(a > PETSC_MUMPS_INT_MAX || a < PETSC_MUMPS_INT_MIN)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"PetscInt too long for PetscMUMPSInt");
 65:   *b = (PetscMUMPSInt)(a);
 66:   return(0);
 67: }

 69: /* Put these utility routines here since they are only used in this file */
 70: PETSC_STATIC_INLINE PetscErrorCode PetscOptionsMUMPSInt_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscMUMPSInt currentvalue,PetscMUMPSInt *value,PetscBool *set,PetscMUMPSInt lb,PetscMUMPSInt ub)
 71: {
 73:   PetscInt       myval;
 74:   PetscBool      myset;
 76:   /* PetscInt's size should be always >= PetscMUMPSInt's. It is safe to call PetscOptionsInt_Private to read a PetscMUMPSInt */
 77:   PetscOptionsInt_Private(PetscOptionsObject,opt,text,man,(PetscInt)currentvalue,&myval,&myset,lb,ub);
 78:   if (myset) {PetscMUMPSIntCast(myval,value);}
 79:   if (set) *set = myset;
 80:   return(0);
 81: }
 82: #define PetscOptionsMUMPSInt(a,b,c,d,e,f) PetscOptionsMUMPSInt_Private(PetscOptionsObject,a,b,c,d,e,f,PETSC_MUMPS_INT_MIN,PETSC_MUMPS_INT_MAX)

 84: /* 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 */
 85: #if defined(PETSC_HAVE_OPENMP_SUPPORT)
 86: #define PetscMUMPS_c(mumps) \
 87:   do { \
 88:     if (mumps->use_petsc_omp_support) { \
 89:       if (mumps->is_omp_master) { \
 90:         PetscOmpCtrlOmpRegionOnMasterBegin(mumps->omp_ctrl); \
 91:         MUMPS_c(&mumps->id); \
 92:         PetscOmpCtrlOmpRegionOnMasterEnd(mumps->omp_ctrl); \
 93:       } \
 94:       PetscOmpCtrlBarrier(mumps->omp_ctrl); \
 95:       /* Global info is same on all processes so we Bcast it within omp_comm. Local info is specific      \
 96:          to processes, so we only Bcast info[1], an error code and leave others (since they do not have   \
 97:          an easy translation between omp_comm and petsc_comm). See MUMPS-5.1.2 manual p82.                   \
 98:          omp_comm is a small shared memory communicator, hence doing multiple Bcast as shown below is OK. \
 99:       */ \
100:       MPI_Bcast(mumps->id.infog, 40,MPIU_MUMPSINT, 0,mumps->omp_comm);  \
101:       MPI_Bcast(mumps->id.rinfog,20,MPIU_REAL,     0,mumps->omp_comm); \
102:       MPI_Bcast(mumps->id.info,  1, MPIU_MUMPSINT, 0,mumps->omp_comm);  \
103:     } else { \
104:       MUMPS_c(&mumps->id); \
105:     } \
106:   } while(0)
107: #else
108: #define PetscMUMPS_c(mumps) \
109:   do { MUMPS_c(&mumps->id); } while (0)
110: #endif

112: /* declare MumpsScalar */
113: #if defined(PETSC_USE_COMPLEX)
114: #if defined(PETSC_USE_REAL_SINGLE)
115: #define MumpsScalar mumps_complex
116: #else
117: #define MumpsScalar mumps_double_complex
118: #endif
119: #else
120: #define MumpsScalar PetscScalar
121: #endif

123: /* macros s.t. indices match MUMPS documentation */
124: #define ICNTL(I) icntl[(I)-1]
125: #define CNTL(I) cntl[(I)-1]
126: #define INFOG(I) infog[(I)-1]
127: #define INFO(I) info[(I)-1]
128: #define RINFOG(I) rinfog[(I)-1]
129: #define RINFO(I) rinfo[(I)-1]

131: typedef struct Mat_MUMPS Mat_MUMPS;
132: struct Mat_MUMPS {
133: #if defined(PETSC_USE_COMPLEX)
134: #if defined(PETSC_USE_REAL_SINGLE)
135:   CMUMPS_STRUC_C id;
136: #else
137:   ZMUMPS_STRUC_C id;
138: #endif
139: #else
140: #if defined(PETSC_USE_REAL_SINGLE)
141:   SMUMPS_STRUC_C id;
142: #else
143:   DMUMPS_STRUC_C id;
144: #endif
145: #endif

147:   MatStructure   matstruc;
148:   PetscMPIInt    myid,petsc_size;
149:   PetscMUMPSInt  *irn,*jcn;             /* the (i,j,v) triplets passed to mumps. */
150:   PetscScalar    *val,*val_alloc;       /* For some matrices, we can directly access their data array without a buffer. For others, we need a buffer. So comes val_alloc. */
151:   PetscInt64     nnz;                   /* number of nonzeros. The type is called selective 64-bit in mumps */
152:   PetscMUMPSInt  sym;
153:   MPI_Comm       mumps_comm;
154:   PetscMUMPSInt  ICNTL9_pre;            /* check if ICNTL(9) is changed from previous MatSolve */
155:   VecScatter     scat_rhs, scat_sol;    /* used by MatSolve() */
156:   Vec            b_seq,x_seq;
157:   PetscInt       ninfo,*info;           /* which INFO to display */
158:   PetscInt       sizeredrhs;
159:   PetscScalar    *schur_sol;
160:   PetscInt       schur_sizesol;
161:   PetscMUMPSInt  *ia_alloc,*ja_alloc;   /* work arrays used for the CSR struct for sparse rhs */
162:   PetscInt64     cur_ilen,cur_jlen;     /* current len of ia_alloc[], ja_alloc[] */
163:   PetscErrorCode (*ConvertToTriples)(Mat,PetscInt,MatReuse,Mat_MUMPS*);

165:   /* stuff used by petsc/mumps OpenMP support*/
166:   PetscBool      use_petsc_omp_support;
167:   PetscOmpCtrl   omp_ctrl;              /* an OpenMP controler that blocked processes will release their CPU (MPI_Barrier does not have this guarantee) */
168:   MPI_Comm       petsc_comm,omp_comm;   /* petsc_comm is petsc matrix's comm */
169:   PetscInt64     *recvcount;            /* a collection of nnz on omp_master */
170:   PetscMPIInt    tag,omp_comm_size;
171:   PetscBool      is_omp_master;         /* is this rank the master of omp_comm */
172:   MPI_Request    *reqs;
173: };

175: /* Cast a 1-based CSR represented by (nrow, ia, ja) of type PetscInt to a CSR of type PetscMUMPSInt.
176:    Here, nrow is number of rows, ia[] is row pointer and ja[] is column indices.
177:  */
178: static PetscErrorCode PetscMUMPSIntCSRCast(Mat_MUMPS *mumps,PetscInt nrow,PetscInt *ia,PetscInt *ja,PetscMUMPSInt **ia_mumps,PetscMUMPSInt **ja_mumps,PetscMUMPSInt *nnz_mumps)
179: {
181:   PetscInt       nnz=ia[nrow]-1; /* mumps uses 1-based indices. Uses PetscInt instead of PetscInt64 since mumps only uses PetscMUMPSInt for rhs */

184: #if defined(PETSC_USE_64BIT_INDICES)
185:   {
186:     PetscInt i;
187:     if (nrow+1 > mumps->cur_ilen) { /* realloc ia_alloc/ja_alloc to fit ia/ja */
188:       PetscFree(mumps->ia_alloc);
189:       PetscMalloc1(nrow+1,&mumps->ia_alloc);
190:       mumps->cur_ilen = nrow+1;
191:     }
192:     if (nnz > mumps->cur_jlen) {
193:       PetscFree(mumps->ja_alloc);
194:       PetscMalloc1(nnz,&mumps->ja_alloc);
195:       mumps->cur_jlen = nnz;
196:     }
197:     for (i=0; i<nrow+1; i++) {PetscMUMPSIntCast(ia[i],&(mumps->ia_alloc[i]));}
198:     for (i=0; i<nnz; i++)    {PetscMUMPSIntCast(ja[i],&(mumps->ja_alloc[i]));}
199:     *ia_mumps = mumps->ia_alloc;
200:     *ja_mumps = mumps->ja_alloc;
201:   }
202: #else
203:   *ia_mumps = ia;
204:   *ja_mumps = ja;
205: #endif
206:   PetscMUMPSIntCast(nnz,nnz_mumps);
207:   return(0);
208: }

210: static PetscErrorCode MatMumpsResetSchur_Private(Mat_MUMPS* mumps)
211: {

215:   PetscFree(mumps->id.listvar_schur);
216:   PetscFree(mumps->id.redrhs);
217:   PetscFree(mumps->schur_sol);
218:   mumps->id.size_schur = 0;
219:   mumps->id.schur_lld  = 0;
220:   mumps->id.ICNTL(19)  = 0;
221:   return(0);
222: }

224: /* solve with rhs in mumps->id.redrhs and return in the same location */
225: static PetscErrorCode MatMumpsSolveSchur_Private(Mat F)
226: {
227:   Mat_MUMPS            *mumps=(Mat_MUMPS*)F->data;
228:   Mat                  S,B,X;
229:   MatFactorSchurStatus schurstatus;
230:   PetscInt             sizesol;
231:   PetscErrorCode       ierr;

234:   MatFactorFactorizeSchurComplement(F);
235:   MatFactorGetSchurComplement(F,&S,&schurstatus);
236:   MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.nrhs,(PetscScalar*)mumps->id.redrhs,&B);
237:   MatSetType(B,((PetscObject)S)->type_name);
238: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
239:   MatBindToCPU(B,S->boundtocpu);
240: #endif
241:   switch (schurstatus) {
242:   case MAT_FACTOR_SCHUR_FACTORED:
243:     MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.nrhs,(PetscScalar*)mumps->id.redrhs,&X);
244:     MatSetType(X,((PetscObject)S)->type_name);
245: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
246:     MatBindToCPU(X,S->boundtocpu);
247: #endif
248:     if (!mumps->id.ICNTL(9)) { /* transpose solve */
249:       MatMatSolveTranspose(S,B,X);
250:     } else {
251:       MatMatSolve(S,B,X);
252:     }
253:     break;
254:   case MAT_FACTOR_SCHUR_INVERTED:
255:     sizesol = mumps->id.nrhs*mumps->id.size_schur;
256:     if (!mumps->schur_sol || sizesol > mumps->schur_sizesol) {
257:       PetscFree(mumps->schur_sol);
258:       PetscMalloc1(sizesol,&mumps->schur_sol);
259:       mumps->schur_sizesol = sizesol;
260:     }
261:     MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.nrhs,mumps->schur_sol,&X);
262:     MatSetType(X,((PetscObject)S)->type_name);
263: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
264:     MatBindToCPU(X,S->boundtocpu);
265: #endif
266:     MatProductCreateWithMat(S,B,NULL,X);
267:     if (!mumps->id.ICNTL(9)) { /* transpose solve */
268:       MatProductSetType(X,MATPRODUCT_AtB);
269:     } else {
270:       MatProductSetType(X,MATPRODUCT_AB);
271:     }
272:     MatProductSetFromOptions(X);
273:     MatProductSymbolic(X);
274:     MatProductNumeric(X);

276:     MatCopy(X,B,SAME_NONZERO_PATTERN);
277:     break;
278:   default:
279:     SETERRQ1(PetscObjectComm((PetscObject)F),PETSC_ERR_SUP,"Unhandled MatFactorSchurStatus %D",F->schur_status);
280:     break;
281:   }
282:   MatFactorRestoreSchurComplement(F,&S,schurstatus);
283:   MatDestroy(&B);
284:   MatDestroy(&X);
285:   return(0);
286: }

288: static PetscErrorCode MatMumpsHandleSchur_Private(Mat F, PetscBool expansion)
289: {
290:   Mat_MUMPS     *mumps=(Mat_MUMPS*)F->data;

294:   if (!mumps->id.ICNTL(19)) { /* do nothing when Schur complement has not been computed */
295:     return(0);
296:   }
297:   if (!expansion) { /* prepare for the condensation step */
298:     PetscInt sizeredrhs = mumps->id.nrhs*mumps->id.size_schur;
299:     /* allocate MUMPS internal array to store reduced right-hand sides */
300:     if (!mumps->id.redrhs || sizeredrhs > mumps->sizeredrhs) {
301:       PetscFree(mumps->id.redrhs);
302:       mumps->id.lredrhs = mumps->id.size_schur;
303:       PetscMalloc1(mumps->id.nrhs*mumps->id.lredrhs,&mumps->id.redrhs);
304:       mumps->sizeredrhs = mumps->id.nrhs*mumps->id.lredrhs;
305:     }
306:     mumps->id.ICNTL(26) = 1; /* condensation phase */
307:   } else { /* prepare for the expansion step */
308:     /* solve Schur complement (this has to be done by the MUMPS user, so basically us) */
309:     MatMumpsSolveSchur_Private(F);
310:     mumps->id.ICNTL(26) = 2; /* expansion phase */
311:     PetscMUMPS_c(mumps);
312:     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));
313:     /* restore defaults */
314:     mumps->id.ICNTL(26) = -1;
315:     /* free MUMPS internal array for redrhs if we have solved for multiple rhs in order to save memory space */
316:     if (mumps->id.nrhs > 1) {
317:       PetscFree(mumps->id.redrhs);
318:       mumps->id.lredrhs = 0;
319:       mumps->sizeredrhs = 0;
320:     }
321:   }
322:   return(0);
323: }

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

328:   input:
329:     A       - matrix in aij,baij or sbaij format
330:     shift   - 0: C style output triple; 1: Fortran style output triple.
331:     reuse   - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple
332:               MAT_REUSE_MATRIX:   only the values in v array are updated
333:   output:
334:     nnz     - dim of r, c, and v (number of local nonzero entries of A)
335:     r, c, v - row and col index, matrix values (matrix triples)

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

341:  */

343: PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,PetscInt shift,MatReuse reuse,Mat_MUMPS *mumps)
344: {
345:   const PetscScalar *av;
346:   const PetscInt    *ai,*aj,*ajj,M=A->rmap->n;
347:   PetscInt64        nz,rnz,i,j,k;
348:   PetscErrorCode    ierr;
349:   PetscMUMPSInt     *row,*col;
350:   Mat_SeqAIJ        *aa=(Mat_SeqAIJ*)A->data;

353:   MatSeqAIJGetArrayRead(A,&av);
354:   mumps->val = (PetscScalar*)av;
355:   if (reuse == MAT_INITIAL_MATRIX) {
356:     nz   = aa->nz;
357:     ai   = aa->i;
358:     aj   = aa->j;
359:     PetscMalloc2(nz,&row,nz,&col);
360:     for (i=k=0; i<M; i++) {
361:       rnz = ai[i+1] - ai[i];
362:       ajj = aj + ai[i];
363:       for (j=0; j<rnz; j++) {
364:         PetscMUMPSIntCast(i+shift,&row[k]);
365:         PetscMUMPSIntCast(ajj[j] + shift,&col[k]);
366:         k++;
367:       }
368:     }
369:     mumps->irn = row;
370:     mumps->jcn = col;
371:     mumps->nnz = nz;
372:   }
373:   MatSeqAIJRestoreArrayRead(A,&av);
374:   return(0);
375: }

377: PetscErrorCode MatConvertToTriples_seqsell_seqaij(Mat A,PetscInt shift,MatReuse reuse,Mat_MUMPS *mumps)
378: {
380:   PetscInt64     nz,i,j,k,r;
381:   Mat_SeqSELL    *a=(Mat_SeqSELL*)A->data;
382:   PetscMUMPSInt  *row,*col;

385:   mumps->val = a->val;
386:   if (reuse == MAT_INITIAL_MATRIX) {
387:     nz   = a->sliidx[a->totalslices];
388:     PetscMalloc2(nz,&row,nz,&col);
389:     for (i=k=0; i<a->totalslices; i++) {
390:       for (j=a->sliidx[i],r=0; j<a->sliidx[i+1]; j++,r=((r+1)&0x07)) {
391:         PetscMUMPSIntCast(8*i+r+shift,&row[k++]);
392:       }
393:     }
394:     for (i=0;i<nz;i++) {PetscMUMPSIntCast(a->colidx[i]+shift,&col[i]);}
395:     mumps->irn = row;
396:     mumps->jcn = col;
397:     mumps->nnz = nz;
398:   }
399:   return(0);
400: }

402: PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,PetscInt shift,MatReuse reuse,Mat_MUMPS *mumps)
403: {
404:   Mat_SeqBAIJ    *aa=(Mat_SeqBAIJ*)A->data;
405:   const PetscInt *ai,*aj,*ajj,bs2 = aa->bs2;
406:   PetscInt64     M,nz,idx=0,rnz,i,j,k,m;
407:   PetscInt       bs;
409:   PetscMUMPSInt  *row,*col;

412:   MatGetBlockSize(A,&bs);
413:   M          = A->rmap->N/bs;
414:   mumps->val = aa->a;
415:   if (reuse == MAT_INITIAL_MATRIX) {
416:     ai   = aa->i; aj = aa->j;
417:     nz   = bs2*aa->nz;
418:     PetscMalloc2(nz,&row,nz,&col);
419:     for (i=0; i<M; i++) {
420:       ajj = aj + ai[i];
421:       rnz = ai[i+1] - ai[i];
422:       for (k=0; k<rnz; k++) {
423:         for (j=0; j<bs; j++) {
424:           for (m=0; m<bs; m++) {
425:             PetscMUMPSIntCast(i*bs + m + shift,&row[idx]);
426:             PetscMUMPSIntCast(bs*ajj[k] + j + shift,&col[idx]);
427:             idx++;
428:           }
429:         }
430:       }
431:     }
432:     mumps->irn = row;
433:     mumps->jcn = col;
434:     mumps->nnz = nz;
435:   }
436:   return(0);
437: }

439: PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,PetscInt shift,MatReuse reuse,Mat_MUMPS *mumps)
440: {
441:   const PetscInt *ai, *aj,*ajj;
442:   PetscInt        bs;
443:   PetscInt64      nz,rnz,i,j,k,m;
444:   PetscErrorCode  ierr;
445:   PetscMUMPSInt   *row,*col;
446:   PetscScalar     *val;
447:   Mat_SeqSBAIJ    *aa=(Mat_SeqSBAIJ*)A->data;
448:   const PetscInt  bs2=aa->bs2,mbs=aa->mbs;
449: #if defined(PETSC_USE_COMPLEX)
450:   PetscBool       hermitian;
451: #endif

454: #if defined(PETSC_USE_COMPLEX)
455:   MatGetOption(A,MAT_HERMITIAN,&hermitian);
456:   if (hermitian) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MUMPS does not support Hermitian symmetric matrices for Choleksy");
457: #endif
458:   ai   = aa->i;
459:   aj   = aa->j;
460:   MatGetBlockSize(A,&bs);
461:   if (reuse == MAT_INITIAL_MATRIX) {
462:     nz   = aa->nz;
463:     PetscMalloc2(bs2*nz,&row,bs2*nz,&col);
464:     if (bs>1) {
465:       PetscMalloc1(bs2*nz,&mumps->val_alloc);
466:       mumps->val = mumps->val_alloc;
467:     } else {
468:       mumps->val = aa->a;
469:     }
470:     mumps->irn = row;
471:     mumps->jcn = col;
472:   } else {
473:     if (bs == 1) mumps->val = aa->a;
474:     row = mumps->irn;
475:     col = mumps->jcn;
476:   }
477:   val = mumps->val;

479:   nz = 0;
480:   if (bs>1) {
481:     for (i=0; i<mbs; i++) {
482:       rnz = ai[i+1] - ai[i];
483:       ajj = aj + ai[i];
484:       for (j=0; j<rnz; j++) {
485:         for (k=0; k<bs; k++) {
486:           for (m=0; m<bs; m++) {
487:             if (ajj[j]>i || k>=m) {
488:               if (reuse == MAT_INITIAL_MATRIX) {
489:                 PetscMUMPSIntCast(i*bs + m + shift,&row[nz]);
490:                 PetscMUMPSIntCast(ajj[j]*bs + k + shift,&col[nz]);
491:               }
492:               val[nz++] = aa->a[(ai[i]+j)*bs2 + m + k*bs];
493:             }
494:           }
495:         }
496:       }
497:     }
498:   } else if (reuse == MAT_INITIAL_MATRIX) {
499:     for (i=0; i<mbs; i++) {
500:       rnz = ai[i+1] - ai[i];
501:       ajj = aj + ai[i];
502:       for (j=0; j<rnz; j++) {
503:         PetscMUMPSIntCast(i+shift,&row[nz]);
504:         PetscMUMPSIntCast(ajj[j] + shift,&col[nz]);
505:         nz++;
506:       }
507:     }
508:     if (nz != aa->nz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Different numbers of nonzeros %D != %D",nz,aa->nz);
509:   }
510:   if (reuse == MAT_INITIAL_MATRIX) mumps->nnz = nz;
511:   return(0);
512: }

514: PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,PetscInt shift,MatReuse reuse,Mat_MUMPS *mumps)
515: {
516:   const PetscInt    *ai,*aj,*ajj,*adiag,M=A->rmap->n;
517:   PetscInt64        nz,rnz,i,j;
518:   const PetscScalar *av,*v1;
519:   PetscScalar       *val;
520:   PetscErrorCode    ierr;
521:   PetscMUMPSInt     *row,*col;
522:   Mat_SeqAIJ        *aa=(Mat_SeqAIJ*)A->data;
523:   PetscBool         missing;
524: #if defined(PETSC_USE_COMPLEX)
525:   PetscBool         hermitian;
526: #endif

529: #if defined(PETSC_USE_COMPLEX)
530:   MatGetOption(A,MAT_HERMITIAN,&hermitian);
531:   if (hermitian) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MUMPS does not support Hermitian symmetric matrices for Choleksy");
532: #endif
533:   MatSeqAIJGetArrayRead(A,&av);
534:   ai    = aa->i; aj = aa->j;
535:   adiag = aa->diag;
536:   MatMissingDiagonal_SeqAIJ(A,&missing,NULL);
537:   if (reuse == MAT_INITIAL_MATRIX) {
538:     /* count nz in the upper triangular part of A */
539:     nz = 0;
540:     if (missing) {
541:       for (i=0; i<M; i++) {
542:         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
543:           for (j=ai[i];j<ai[i+1];j++) {
544:             if (aj[j] < i) continue;
545:             nz++;
546:           }
547:         } else {
548:           nz += ai[i+1] - adiag[i];
549:         }
550:       }
551:     } else {
552:       for (i=0; i<M; i++) nz += ai[i+1] - adiag[i];
553:     }
554:     PetscMalloc2(nz,&row,nz,&col);
555:     PetscMalloc1(nz,&val);
556:     mumps->nnz = nz;
557:     mumps->irn = row;
558:     mumps->jcn = col;
559:     mumps->val = mumps->val_alloc = val;

561:     nz = 0;
562:     if (missing) {
563:       for (i=0; i<M; i++) {
564:         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
565:           for (j=ai[i];j<ai[i+1];j++) {
566:             if (aj[j] < i) continue;
567:             PetscMUMPSIntCast(i+shift,&row[nz]);
568:             PetscMUMPSIntCast(aj[j]+shift,&col[nz]);
569:             val[nz] = av[j];
570:             nz++;
571:           }
572:         } else {
573:           rnz = ai[i+1] - adiag[i];
574:           ajj = aj + adiag[i];
575:           v1  = av + adiag[i];
576:           for (j=0; j<rnz; j++) {
577:             PetscMUMPSIntCast(i+shift,&row[nz]);
578:             PetscMUMPSIntCast(ajj[j] + shift,&col[nz]);
579:             val[nz++] = v1[j];
580:           }
581:         }
582:       }
583:     } else {
584:       for (i=0; i<M; i++) {
585:         rnz = ai[i+1] - adiag[i];
586:         ajj = aj + adiag[i];
587:         v1  = av + adiag[i];
588:         for (j=0; j<rnz; j++) {
589:           PetscMUMPSIntCast(i+shift,&row[nz]);
590:           PetscMUMPSIntCast(ajj[j] + shift,&col[nz]);
591:           val[nz++] = v1[j];
592:         }
593:       }
594:     }
595:   } else {
596:     nz = 0;
597:     val = mumps->val;
598:     if (missing) {
599:       for (i=0; i <M; i++) {
600:         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
601:           for (j=ai[i];j<ai[i+1];j++) {
602:             if (aj[j] < i) continue;
603:             val[nz++] = av[j];
604:           }
605:         } else {
606:           rnz = ai[i+1] - adiag[i];
607:           v1  = av + adiag[i];
608:           for (j=0; j<rnz; j++) {
609:             val[nz++] = v1[j];
610:           }
611:         }
612:       }
613:     } else {
614:       for (i=0; i <M; i++) {
615:         rnz = ai[i+1] - adiag[i];
616:         v1  = av + adiag[i];
617:         for (j=0; j<rnz; j++) {
618:           val[nz++] = v1[j];
619:         }
620:       }
621:     }
622:   }
623:   MatSeqAIJRestoreArrayRead(A,&av);
624:   return(0);
625: }

627: PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,PetscInt shift,MatReuse reuse,Mat_MUMPS *mumps)
628: {
629:   PetscErrorCode    ierr;
630:   const PetscInt    *ai,*aj,*bi,*bj,*garray,*ajj,*bjj;
631:   PetscInt          bs;
632:   PetscInt64        rstart,nz,i,j,k,m,jj,irow,countA,countB;
633:   PetscMUMPSInt     *row,*col;
634:   const PetscScalar *av,*bv,*v1,*v2;
635:   PetscScalar       *val;
636:   Mat_MPISBAIJ      *mat = (Mat_MPISBAIJ*)A->data;
637:   Mat_SeqSBAIJ      *aa  = (Mat_SeqSBAIJ*)(mat->A)->data;
638:   Mat_SeqBAIJ       *bb  = (Mat_SeqBAIJ*)(mat->B)->data;
639:   const PetscInt    bs2=aa->bs2,mbs=aa->mbs;
640: #if defined(PETSC_USE_COMPLEX)
641:   PetscBool         hermitian;
642: #endif

645: #if defined(PETSC_USE_COMPLEX)
646:   MatGetOption(A,MAT_HERMITIAN,&hermitian);
647:   if (hermitian) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MUMPS does not support Hermitian symmetric matrices for Choleksy");
648: #endif
649:   MatGetBlockSize(A,&bs);
650:   rstart = A->rmap->rstart;
651:   ai = aa->i;
652:   aj = aa->j;
653:   bi = bb->i;
654:   bj = bb->j;
655:   av = aa->a;
656:   bv = bb->a;

658:   garray = mat->garray;

660:   if (reuse == MAT_INITIAL_MATRIX) {
661:     nz   = (aa->nz+bb->nz)*bs2; /* just a conservative estimate */
662:     PetscMalloc2(nz,&row,nz,&col);
663:     PetscMalloc1(nz,&val);
664:     /* can not decide the exact mumps->nnz now because of the SBAIJ */
665:     mumps->irn = row;
666:     mumps->jcn = col;
667:     mumps->val = mumps->val_alloc = val;
668:   } else {
669:     val = mumps->val;
670:   }

672:   jj = 0; irow = rstart;
673:   for (i=0; i<mbs; i++) {
674:     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
675:     countA = ai[i+1] - ai[i];
676:     countB = bi[i+1] - bi[i];
677:     bjj    = bj + bi[i];
678:     v1     = av + ai[i]*bs2;
679:     v2     = bv + bi[i]*bs2;

681:     if (bs>1) {
682:       /* A-part */
683:       for (j=0; j<countA; j++) {
684:         for (k=0; k<bs; k++) {
685:           for (m=0; m<bs; m++) {
686:             if (rstart + ajj[j]*bs>irow || k>=m) {
687:               if (reuse == MAT_INITIAL_MATRIX) {
688:                 PetscMUMPSIntCast(irow + m + shift,&row[jj]);
689:                 PetscMUMPSIntCast(rstart + ajj[j]*bs + k + shift,&col[jj]);
690:               }
691:               val[jj++] = v1[j*bs2 + m + k*bs];
692:             }
693:           }
694:         }
695:       }

697:       /* B-part */
698:       for (j=0; j < countB; j++) {
699:         for (k=0; k<bs; k++) {
700:           for (m=0; m<bs; m++) {
701:             if (reuse == MAT_INITIAL_MATRIX) {
702:               PetscMUMPSIntCast(irow + m + shift,&row[jj]);
703:               PetscMUMPSIntCast(garray[bjj[j]]*bs + k + shift,&col[jj]);
704:             }
705:             val[jj++] = v2[j*bs2 + m + k*bs];
706:           }
707:         }
708:       }
709:     } else {
710:       /* A-part */
711:       for (j=0; j<countA; j++) {
712:         if (reuse == MAT_INITIAL_MATRIX) {
713:           PetscMUMPSIntCast(irow + shift,&row[jj]);
714:           PetscMUMPSIntCast(rstart + ajj[j] + shift,&col[jj]);
715:         }
716:         val[jj++] = v1[j];
717:       }

719:       /* B-part */
720:       for (j=0; j < countB; j++) {
721:         if (reuse == MAT_INITIAL_MATRIX) {
722:           PetscMUMPSIntCast(irow + shift,&row[jj]);
723:           PetscMUMPSIntCast(garray[bjj[j]] + shift,&col[jj]);
724:         }
725:         val[jj++] = v2[j];
726:       }
727:     }
728:     irow+=bs;
729:   }
730:   mumps->nnz = jj;
731:   return(0);
732: }

734: PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,PetscInt shift,MatReuse reuse,Mat_MUMPS *mumps)
735: {
736:   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
737:   PetscErrorCode    ierr;
738:   PetscInt64        rstart,nz,i,j,jj,irow,countA,countB;
739:   PetscMUMPSInt     *row,*col;
740:   const PetscScalar *av, *bv,*v1,*v2;
741:   PetscScalar       *val;
742:   Mat               Ad,Ao;
743:   Mat_SeqAIJ        *aa;
744:   Mat_SeqAIJ        *bb;

747:   MatMPIAIJGetSeqAIJ(A,&Ad,&Ao,&garray);
748:   MatSeqAIJGetArrayRead(Ad,&av);
749:   MatSeqAIJGetArrayRead(Ao,&bv);

751:   aa = (Mat_SeqAIJ*)(Ad)->data;
752:   bb = (Mat_SeqAIJ*)(Ao)->data;
753:   ai = aa->i;
754:   aj = aa->j;
755:   bi = bb->i;
756:   bj = bb->j;

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

760:   if (reuse == MAT_INITIAL_MATRIX) {
761:     nz   = (PetscInt64)aa->nz + bb->nz; /* make sure the sum won't overflow PetscInt */
762:     PetscMalloc2(nz,&row,nz,&col);
763:     PetscMalloc1(nz,&val);
764:     mumps->nnz = nz;
765:     mumps->irn = row;
766:     mumps->jcn = col;
767:     mumps->val = mumps->val_alloc = val;
768:   } else {
769:     val = mumps->val;
770:   }

772:   jj = 0; irow = rstart;
773:   for (i=0; i<m; i++) {
774:     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
775:     countA = ai[i+1] - ai[i];
776:     countB = bi[i+1] - bi[i];
777:     bjj    = bj + bi[i];
778:     v1     = av + ai[i];
779:     v2     = bv + bi[i];

781:     /* A-part */
782:     for (j=0; j<countA; j++) {
783:       if (reuse == MAT_INITIAL_MATRIX) {
784:         PetscMUMPSIntCast(irow + shift,&row[jj]);
785:         PetscMUMPSIntCast(rstart + ajj[j] + shift,&col[jj]);
786:       }
787:       val[jj++] = v1[j];
788:     }

790:     /* B-part */
791:     for (j=0; j < countB; j++) {
792:       if (reuse == MAT_INITIAL_MATRIX) {
793:         PetscMUMPSIntCast(irow + shift,&row[jj]);
794:         PetscMUMPSIntCast(garray[bjj[j]] + shift,&col[jj]);
795:       }
796:       val[jj++] = v2[j];
797:     }
798:     irow++;
799:   }
800:   MatSeqAIJRestoreArrayRead(Ad,&av);
801:   MatSeqAIJRestoreArrayRead(Ao,&bv);
802:   return(0);
803: }

805: PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,PetscInt shift,MatReuse reuse,Mat_MUMPS *mumps)
806: {
807:   Mat_MPIBAIJ       *mat    = (Mat_MPIBAIJ*)A->data;
808:   Mat_SeqBAIJ       *aa     = (Mat_SeqBAIJ*)(mat->A)->data;
809:   Mat_SeqBAIJ       *bb     = (Mat_SeqBAIJ*)(mat->B)->data;
810:   const PetscInt    *ai     = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj;
811:   const PetscInt    *garray = mat->garray,mbs=mat->mbs,rstart=A->rmap->rstart;
812:   const PetscInt    bs2=mat->bs2;
813:   PetscErrorCode    ierr;
814:   PetscInt          bs;
815:   PetscInt64        nz,i,j,k,n,jj,irow,countA,countB,idx;
816:   PetscMUMPSInt     *row,*col;
817:   const PetscScalar *av=aa->a, *bv=bb->a,*v1,*v2;
818:   PetscScalar       *val;

821:   MatGetBlockSize(A,&bs);
822:   if (reuse == MAT_INITIAL_MATRIX) {
823:     nz   = bs2*(aa->nz + bb->nz);
824:     PetscMalloc2(nz,&row,nz,&col);
825:     PetscMalloc1(nz,&val);
826:     mumps->nnz = nz;
827:     mumps->irn = row;
828:     mumps->jcn = col;
829:     mumps->val = mumps->val_alloc = val;
830:   } else {
831:     val = mumps->val;
832:   }

834:   jj = 0; irow = rstart;
835:   for (i=0; i<mbs; i++) {
836:     countA = ai[i+1] - ai[i];
837:     countB = bi[i+1] - bi[i];
838:     ajj    = aj + ai[i];
839:     bjj    = bj + bi[i];
840:     v1     = av + bs2*ai[i];
841:     v2     = bv + bs2*bi[i];

843:     idx = 0;
844:     /* A-part */
845:     for (k=0; k<countA; k++) {
846:       for (j=0; j<bs; j++) {
847:         for (n=0; n<bs; n++) {
848:           if (reuse == MAT_INITIAL_MATRIX) {
849:             PetscMUMPSIntCast(irow + n + shift,&row[jj]);
850:             PetscMUMPSIntCast(rstart + bs*ajj[k] + j + shift,&col[jj]);
851:           }
852:           val[jj++] = v1[idx++];
853:         }
854:       }
855:     }

857:     idx = 0;
858:     /* B-part */
859:     for (k=0; k<countB; k++) {
860:       for (j=0; j<bs; j++) {
861:         for (n=0; n<bs; n++) {
862:           if (reuse == MAT_INITIAL_MATRIX) {
863:             PetscMUMPSIntCast(irow + n + shift,&row[jj]);
864:             PetscMUMPSIntCast(bs*garray[bjj[k]] + j + shift,&col[jj]);
865:           }
866:           val[jj++] = v2[idx++];
867:         }
868:       }
869:     }
870:     irow += bs;
871:   }
872:   return(0);
873: }

875: PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,PetscInt shift,MatReuse reuse,Mat_MUMPS *mumps)
876: {
877:   const PetscInt    *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
878:   PetscErrorCode    ierr;
879:   PetscInt64        rstart,nz,nza,nzb,i,j,jj,irow,countA,countB;
880:   PetscMUMPSInt     *row,*col;
881:   const PetscScalar *av, *bv,*v1,*v2;
882:   PetscScalar       *val;
883:   Mat               Ad,Ao;
884:   Mat_SeqAIJ        *aa;
885:   Mat_SeqAIJ        *bb;
886: #if defined(PETSC_USE_COMPLEX)
887:   PetscBool         hermitian;
888: #endif

891: #if defined(PETSC_USE_COMPLEX)
892:   MatGetOption(A,MAT_HERMITIAN,&hermitian);
893:   if (hermitian) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MUMPS does not support Hermitian symmetric matrices for Choleksy");
894: #endif
895:   MatMPIAIJGetSeqAIJ(A,&Ad,&Ao,&garray);
896:   MatSeqAIJGetArrayRead(Ad,&av);
897:   MatSeqAIJGetArrayRead(Ao,&bv);

899:   aa    = (Mat_SeqAIJ*)(Ad)->data;
900:   bb    = (Mat_SeqAIJ*)(Ao)->data;
901:   ai    = aa->i;
902:   aj    = aa->j;
903:   adiag = aa->diag;
904:   bi    = bb->i;
905:   bj    = bb->j;

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

909:   if (reuse == MAT_INITIAL_MATRIX) {
910:     nza = 0;    /* num of upper triangular entries in mat->A, including diagonals */
911:     nzb = 0;    /* num of upper triangular entries in mat->B */
912:     for (i=0; i<m; i++) {
913:       nza   += (ai[i+1] - adiag[i]);
914:       countB = bi[i+1] - bi[i];
915:       bjj    = bj + bi[i];
916:       for (j=0; j<countB; j++) {
917:         if (garray[bjj[j]] > rstart) nzb++;
918:       }
919:     }

921:     nz   = nza + nzb; /* total nz of upper triangular part of mat */
922:     PetscMalloc2(nz,&row,nz,&col);
923:     PetscMalloc1(nz,&val);
924:     mumps->nnz = nz;
925:     mumps->irn = row;
926:     mumps->jcn = col;
927:     mumps->val = mumps->val_alloc = val;
928:   } else {
929:     val = mumps->val;
930:   }

932:   jj = 0; irow = rstart;
933:   for (i=0; i<m; i++) {
934:     ajj    = aj + adiag[i];                 /* ptr to the beginning of the diagonal of this row */
935:     v1     = av + adiag[i];
936:     countA = ai[i+1] - adiag[i];
937:     countB = bi[i+1] - bi[i];
938:     bjj    = bj + bi[i];
939:     v2     = bv + bi[i];

941:     /* A-part */
942:     for (j=0; j<countA; j++) {
943:       if (reuse == MAT_INITIAL_MATRIX) {
944:         PetscMUMPSIntCast(irow + shift,&row[jj]);
945:         PetscMUMPSIntCast(rstart + ajj[j] + shift,&col[jj]);
946:       }
947:       val[jj++] = v1[j];
948:     }

950:     /* B-part */
951:     for (j=0; j < countB; j++) {
952:       if (garray[bjj[j]] > rstart) {
953:         if (reuse == MAT_INITIAL_MATRIX) {
954:           PetscMUMPSIntCast(irow + shift,&row[jj]);
955:           PetscMUMPSIntCast(garray[bjj[j]] + shift,&col[jj]);
956:         }
957:         val[jj++] = v2[j];
958:       }
959:     }
960:     irow++;
961:   }
962:   MatSeqAIJRestoreArrayRead(Ad,&av);
963:   MatSeqAIJRestoreArrayRead(Ao,&bv);
964:   return(0);
965: }

967: PetscErrorCode MatDestroy_MUMPS(Mat A)
968: {
970:   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->data;

973:   PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);
974:   VecScatterDestroy(&mumps->scat_rhs);
975:   VecScatterDestroy(&mumps->scat_sol);
976:   VecDestroy(&mumps->b_seq);
977:   VecDestroy(&mumps->x_seq);
978:   PetscFree(mumps->id.perm_in);
979:   PetscFree2(mumps->irn,mumps->jcn);
980:   PetscFree(mumps->val_alloc);
981:   PetscFree(mumps->info);
982:   MatMumpsResetSchur_Private(mumps);
983:   mumps->id.job = JOB_END;
984:   PetscMUMPS_c(mumps);
985:   if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in MatDestroy_MUMPS: INFOG(1)=%d\n",mumps->id.INFOG(1));
986: #if defined(PETSC_HAVE_OPENMP_SUPPORT)
987:   if (mumps->use_petsc_omp_support) { PetscOmpCtrlDestroy(&mumps->omp_ctrl); }
988: #endif
989:   PetscFree(mumps->ia_alloc);
990:   PetscFree(mumps->ja_alloc);
991:   PetscFree(mumps->recvcount);
992:   PetscFree(mumps->reqs);
993:   PetscFree(A->data);

995:   /* clear composed functions */
996:   PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverType_C",NULL);
997:   PetscObjectComposeFunction((PetscObject)A,"MatFactorSetSchurIS_C",NULL);
998:   PetscObjectComposeFunction((PetscObject)A,"MatFactorCreateSchurComplement_C",NULL);
999:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetIcntl_C",NULL);
1000:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetIcntl_C",NULL);
1001:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetCntl_C",NULL);
1002:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetCntl_C",NULL);
1003:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfo_C",NULL);
1004:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfog_C",NULL);
1005:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfo_C",NULL);
1006:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfog_C",NULL);
1007:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInverse_C",NULL);
1008:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInverseTranspose_C",NULL);
1009:   return(0);
1010: }

1012: PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x)
1013: {
1014:   Mat_MUMPS        *mumps=(Mat_MUMPS*)A->data;
1015:   PetscScalar      *array;
1016:   Vec              b_seq;
1017:   IS               is_iden,is_petsc;
1018:   PetscErrorCode   ierr;
1019:   PetscInt         i;
1020:   PetscBool        second_solve = PETSC_FALSE;
1021:   static PetscBool cite1 = PETSC_FALSE,cite2 = PETSC_FALSE;

1024:   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);
1025:   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);

1027:   if (A->factorerrortype) {
1028:     PetscInfo2(A,"MatSolve is called with singular matrix factor, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1029:     VecSetInf(x);
1030:     return(0);
1031:   }

1033:   mumps->id.ICNTL(20) = 0; /* dense RHS */
1034:   mumps->id.nrhs      = 1;
1035:   b_seq               = mumps->b_seq;
1036:   if (mumps->petsc_size > 1) {
1037:     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
1038:     VecScatterBegin(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);
1039:     VecScatterEnd(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);
1040:     if (!mumps->myid) {VecGetArray(b_seq,&array);}
1041:   } else {  /* petsc_size == 1 */
1042:     VecCopy(b,x);
1043:     VecGetArray(x,&array);
1044:   }
1045:   if (!mumps->myid) { /* define rhs on the host */
1046:     mumps->id.nrhs = 1;
1047:     mumps->id.rhs = (MumpsScalar*)array;
1048:   }

1050:   /*
1051:      handle condensation step of Schur complement (if any)
1052:      We set by default ICNTL(26) == -1 when Schur indices have been provided by the user.
1053:      According to MUMPS (5.0.0) manual, any value should be harmful during the factorization phase
1054:      Unless the user provides a valid value for ICNTL(26), MatSolve and MatMatSolve routines solve the full system.
1055:      This requires an extra call to PetscMUMPS_c and the computation of the factors for S
1056:   */
1057:   if (mumps->id.size_schur > 0 && (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2)) {
1058:     if (mumps->petsc_size > 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Parallel Schur complements not yet supported from PETSc\n");
1059:     second_solve = PETSC_TRUE;
1060:     MatMumpsHandleSchur_Private(A,PETSC_FALSE);
1061:   }
1062:   /* solve phase */
1063:   /*-------------*/
1064:   mumps->id.job = JOB_SOLVE;
1065:   PetscMUMPS_c(mumps);
1066:   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));

1068:   /* handle expansion step of Schur complement (if any) */
1069:   if (second_solve) {
1070:     MatMumpsHandleSchur_Private(A,PETSC_TRUE);
1071:   }

1073:   if (mumps->petsc_size > 1) { /* convert mumps distributed solution to petsc mpi x */
1074:     if (mumps->scat_sol && mumps->ICNTL9_pre != mumps->id.ICNTL(9)) {
1075:       /* when id.ICNTL(9) changes, the contents of lsol_loc may change (not its size, lsol_loc), recreates scat_sol */
1076:       VecScatterDestroy(&mumps->scat_sol);
1077:     }
1078:     if (!mumps->scat_sol) { /* create scatter scat_sol */
1079:       PetscInt *isol2_loc=NULL;
1080:       ISCreateStride(PETSC_COMM_SELF,mumps->id.lsol_loc,0,1,&is_iden); /* from */
1081:       PetscMalloc1(mumps->id.lsol_loc,&isol2_loc);
1082:       for (i=0; i<mumps->id.lsol_loc; i++) isol2_loc[i] = mumps->id.isol_loc[i]-1; /* change Fortran style to C style */
1083:       ISCreateGeneral(PETSC_COMM_SELF,mumps->id.lsol_loc,isol2_loc,PETSC_OWN_POINTER,&is_petsc);  /* to */
1084:       VecScatterCreate(mumps->x_seq,is_iden,x,is_petsc,&mumps->scat_sol);
1085:       ISDestroy(&is_iden);
1086:       ISDestroy(&is_petsc);
1087:       mumps->ICNTL9_pre = mumps->id.ICNTL(9); /* save current value of id.ICNTL(9) */
1088:     }

1090:     VecScatterBegin(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);
1091:     VecScatterEnd(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);
1092:   }

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

1097:   PetscLogFlops(2.0*mumps->id.RINFO(3));
1098:   return(0);
1099: }

1101: PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x)
1102: {
1103:   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->data;

1107:   mumps->id.ICNTL(9) = 0;
1108:   MatSolve_MUMPS(A,b,x);
1109:   mumps->id.ICNTL(9) = 1;
1110:   return(0);
1111: }

1113: PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X)
1114: {
1115:   PetscErrorCode    ierr;
1116:   Mat               Bt = NULL;
1117:   PetscBool         denseX,denseB,flg,flgT;
1118:   Mat_MUMPS         *mumps=(Mat_MUMPS*)A->data;
1119:   PetscInt          i,nrhs,M;
1120:   PetscScalar       *array;
1121:   const PetscScalar *rbray;
1122:   PetscInt          lsol_loc,nlsol_loc,*idxx,iidx = 0;
1123:   PetscMUMPSInt     *isol_loc,*isol_loc_save;
1124:   PetscScalar       *bray,*sol_loc,*sol_loc_save;
1125:   IS                is_to,is_from;
1126:   PetscInt          k,proc,j,m,myrstart;
1127:   const PetscInt    *rstart;
1128:   Vec               v_mpi,b_seq,msol_loc;
1129:   VecScatter        scat_rhs,scat_sol;
1130:   PetscScalar       *aa;
1131:   PetscInt          spnr,*ia,*ja;
1132:   Mat_MPIAIJ        *b = NULL;

1135:   PetscObjectTypeCompareAny((PetscObject)X,&denseX,MATSEQDENSE,MATMPIDENSE,NULL);
1136:   if (!denseX) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");

1138:   PetscObjectTypeCompareAny((PetscObject)B,&denseB,MATSEQDENSE,MATMPIDENSE,NULL);
1139:   if (denseB) {
1140:     if (B->rmap->n != X->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix B and X must have same row distribution");
1141:     mumps->id.ICNTL(20)= 0; /* dense RHS */
1142:   } else { /* sparse B */
1143:     if (X == B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_IDN,"X and B must be different matrices");
1144:     PetscObjectTypeCompare((PetscObject)B,MATTRANSPOSEMAT,&flgT);
1145:     if (flgT) { /* input B is transpose of actural RHS matrix,
1146:                  because mumps requires sparse compressed COLUMN storage! See MatMatTransposeSolve_MUMPS() */
1147:       MatTransposeGetMat(B,&Bt);
1148:     } else SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B must be MATTRANSPOSEMAT matrix");
1149:     mumps->id.ICNTL(20)= 1; /* sparse RHS */
1150:   }

1152:   MatGetSize(B,&M,&nrhs);
1153:   mumps->id.nrhs = nrhs;
1154:   mumps->id.lrhs = M;
1155:   mumps->id.rhs  = NULL;

1157:   if (mumps->petsc_size == 1) {
1158:     PetscScalar *aa;
1159:     PetscInt    spnr,*ia,*ja;
1160:     PetscBool   second_solve = PETSC_FALSE;

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

1165:     if (denseB) {
1166:       /* copy B to X */
1167:       MatDenseGetArrayRead(B,&rbray);
1168:       PetscArraycpy(array,rbray,M*nrhs);
1169:       MatDenseRestoreArrayRead(B,&rbray);
1170:     } else { /* sparse B */
1171:       MatSeqAIJGetArray(Bt,&aa);
1172:       MatGetRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&flg);
1173:       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot get IJ structure");
1174:       PetscMUMPSIntCSRCast(mumps,spnr,ia,ja,&mumps->id.irhs_ptr,&mumps->id.irhs_sparse,&mumps->id.nz_rhs);
1175:       mumps->id.rhs_sparse  = (MumpsScalar*)aa;
1176:     }
1177:     /* handle condensation step of Schur complement (if any) */
1178:     if (mumps->id.size_schur > 0 && (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2)) {
1179:       second_solve = PETSC_TRUE;
1180:       MatMumpsHandleSchur_Private(A,PETSC_FALSE);
1181:     }
1182:     /* solve phase */
1183:     /*-------------*/
1184:     mumps->id.job = JOB_SOLVE;
1185:     PetscMUMPS_c(mumps);
1186:     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));

1188:     /* handle expansion step of Schur complement (if any) */
1189:     if (second_solve) {
1190:       MatMumpsHandleSchur_Private(A,PETSC_TRUE);
1191:     }
1192:     if (!denseB) { /* sparse B */
1193:       MatSeqAIJRestoreArray(Bt,&aa);
1194:       MatRestoreRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&flg);
1195:       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot restore IJ structure");
1196:     }
1197:     MatDenseRestoreArray(X,&array);
1198:     return(0);
1199:   }

1201:   /*--------- parallel case: MUMPS requires rhs B to be centralized on the host! --------*/
1202:   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");

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

1208:   lsol_loc  = mumps->id.lsol_loc;
1209:   nlsol_loc = nrhs*lsol_loc;     /* length of sol_loc */
1210:   PetscMalloc2(nlsol_loc,&sol_loc,lsol_loc,&isol_loc);
1211:   mumps->id.sol_loc  = (MumpsScalar*)sol_loc;
1212:   mumps->id.isol_loc = isol_loc;

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

1216:   if (denseB) { /* dense B */
1217:     /* TODO: Because of non-contiguous indices, the created vecscatter scat_rhs is not done in MPI_Gather, resulting in
1218:        very inefficient communication. An optimization is to use VecScatterCreateToZero to gather B to rank 0. Then on rank
1219:        0, re-arrange B into desired order, which is a local operation.
1220:      */

1222:     /* scatter v_mpi to b_seq because MUMPS only supports centralized rhs */
1223:     /* wrap dense rhs matrix B into a vector v_mpi */
1224:     MatGetLocalSize(B,&m,NULL);
1225:     MatDenseGetArray(B,&bray);
1226:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)B),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi);
1227:     MatDenseRestoreArray(B,&bray);

1229:     /* scatter v_mpi to b_seq in proc[0]. MUMPS requires rhs to be centralized on the host! */
1230:     if (!mumps->myid) {
1231:       PetscInt *idx;
1232:       /* idx: maps from k-th index of v_mpi to (i,j)-th global entry of B */
1233:       PetscMalloc1(nrhs*M,&idx);
1234:       MatGetOwnershipRanges(B,&rstart);
1235:       k = 0;
1236:       for (proc=0; proc<mumps->petsc_size; proc++){
1237:         for (j=0; j<nrhs; j++){
1238:           for (i=rstart[proc]; i<rstart[proc+1]; i++) idx[k++] = j*M + i;
1239:         }
1240:       }

1242:       VecCreateSeq(PETSC_COMM_SELF,nrhs*M,&b_seq);
1243:       ISCreateGeneral(PETSC_COMM_SELF,nrhs*M,idx,PETSC_OWN_POINTER,&is_to);
1244:       ISCreateStride(PETSC_COMM_SELF,nrhs*M,0,1,&is_from);
1245:     } else {
1246:       VecCreateSeq(PETSC_COMM_SELF,0,&b_seq);
1247:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_to);
1248:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_from);
1249:     }
1250:     VecScatterCreate(v_mpi,is_from,b_seq,is_to,&scat_rhs);
1251:     VecScatterBegin(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);
1252:     ISDestroy(&is_to);
1253:     ISDestroy(&is_from);
1254:     VecScatterEnd(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);

1256:     if (!mumps->myid) { /* define rhs on the host */
1257:       VecGetArray(b_seq,&bray);
1258:       mumps->id.rhs = (MumpsScalar*)bray;
1259:       VecRestoreArray(b_seq,&bray);
1260:     }

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

1265:     /* wrap dense X into a vector v_mpi */
1266:     MatGetLocalSize(X,&m,NULL);
1267:     MatDenseGetArray(X,&bray);
1268:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)X),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi);
1269:     MatDenseRestoreArray(X,&bray);

1271:     if (!mumps->myid) {
1272:       MatSeqAIJGetArray(b->A,&aa);
1273:       MatGetRowIJ(b->A,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&flg);
1274:       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot get IJ structure");
1275:       PetscMUMPSIntCSRCast(mumps,spnr,ia,ja,&mumps->id.irhs_ptr,&mumps->id.irhs_sparse,&mumps->id.nz_rhs);
1276:       mumps->id.rhs_sparse  = (MumpsScalar*)aa;
1277:     } else {
1278:       mumps->id.irhs_ptr    = NULL;
1279:       mumps->id.irhs_sparse = NULL;
1280:       mumps->id.nz_rhs      = 0;
1281:       mumps->id.rhs_sparse  = NULL;
1282:     }
1283:   }

1285:   /* solve phase */
1286:   /*-------------*/
1287:   mumps->id.job = JOB_SOLVE;
1288:   PetscMUMPS_c(mumps);
1289:   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));

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

1295:   /* create scatter scat_sol */
1296:   MatGetOwnershipRanges(X,&rstart);
1297:   /* iidx: index for scatter mumps solution to petsc X */

1299:   ISCreateStride(PETSC_COMM_SELF,nlsol_loc,0,1,&is_from);
1300:   PetscMalloc1(nlsol_loc,&idxx);
1301:   for (i=0; i<lsol_loc; i++) {
1302:     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 */

1304:     for (proc=0; proc<mumps->petsc_size; proc++){
1305:       if (isol_loc[i] >= rstart[proc] && isol_loc[i] < rstart[proc+1]) {
1306:         myrstart = rstart[proc];
1307:         k        = isol_loc[i] - myrstart;        /* local index on 1st column of petsc vector X */
1308:         iidx     = k + myrstart*nrhs;             /* maps mumps isol_loc[i] to petsc index in X */
1309:         m        = rstart[proc+1] - rstart[proc]; /* rows of X for this proc */
1310:         break;
1311:       }
1312:     }

1314:     for (j=0; j<nrhs; j++) idxx[i+j*lsol_loc] = iidx + j*m;
1315:   }
1316:   ISCreateGeneral(PETSC_COMM_SELF,nlsol_loc,idxx,PETSC_COPY_VALUES,&is_to);
1317:   VecScatterCreate(msol_loc,is_from,v_mpi,is_to,&scat_sol);
1318:   VecScatterBegin(scat_sol,msol_loc,v_mpi,INSERT_VALUES,SCATTER_FORWARD);
1319:   ISDestroy(&is_from);
1320:   ISDestroy(&is_to);
1321:   VecScatterEnd(scat_sol,msol_loc,v_mpi,INSERT_VALUES,SCATTER_FORWARD);
1322:   MatDenseRestoreArray(X,&array);

1324:   /* free spaces */
1325:   mumps->id.sol_loc  = (MumpsScalar*)sol_loc_save;
1326:   mumps->id.isol_loc = isol_loc_save;

1328:   PetscFree2(sol_loc,isol_loc);
1329:   PetscFree(idxx);
1330:   VecDestroy(&msol_loc);
1331:   VecDestroy(&v_mpi);
1332:   if (!denseB) {
1333:     if (!mumps->myid) {
1334:       b = (Mat_MPIAIJ*)Bt->data;
1335:       MatSeqAIJRestoreArray(b->A,&aa);
1336:       MatRestoreRowIJ(b->A,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&flg);
1337:       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot restore IJ structure");
1338:     }
1339:   } else {
1340:     VecDestroy(&b_seq);
1341:     VecScatterDestroy(&scat_rhs);
1342:   }
1343:   VecScatterDestroy(&scat_sol);
1344:   PetscLogFlops(2.0*nrhs*mumps->id.RINFO(3));
1345:   return(0);
1346: }

1348: PetscErrorCode MatMatTransposeSolve_MUMPS(Mat A,Mat Bt,Mat X)
1349: {
1351:   PetscBool      flg;
1352:   Mat            B;

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

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

1361:   MatMatSolve_MUMPS(A,B,X);
1362:   MatDestroy(&B);
1363:   return(0);
1364: }

1366: #if !defined(PETSC_USE_COMPLEX)
1367: /*
1368:   input:
1369:    F:        numeric factor
1370:   output:
1371:    nneg:     total number of negative pivots
1372:    nzero:    total number of zero pivots
1373:    npos:     (global dimension of F) - nneg - nzero
1374: */
1375: PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,PetscInt *nneg,PetscInt *nzero,PetscInt *npos)
1376: {
1377:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1379:   PetscMPIInt    size;

1382:   MPI_Comm_size(PetscObjectComm((PetscObject)F),&size);
1383:   /* 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 */
1384:   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));

1386:   if (nneg) *nneg = mumps->id.INFOG(12);
1387:   if (nzero || npos) {
1388:     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");
1389:     if (nzero) *nzero = mumps->id.INFOG(28);
1390:     if (npos) *npos   = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28));
1391:   }
1392:   return(0);
1393: }
1394: #endif

1396: PetscErrorCode MatMumpsGatherNonzerosOnMaster(MatReuse reuse,Mat_MUMPS *mumps)
1397: {
1399:   PetscInt       i,nreqs;
1400:   PetscMUMPSInt  *irn,*jcn;
1401:   PetscMPIInt    count;
1402:   PetscInt64     totnnz,remain;
1403:   const PetscInt osize=mumps->omp_comm_size;
1404:   PetscScalar    *val;

1407:   if (osize > 1) {
1408:     if (reuse == MAT_INITIAL_MATRIX) {
1409:       /* master first gathers counts of nonzeros to receive */
1410:       if (mumps->is_omp_master) {PetscMalloc1(osize,&mumps->recvcount);}
1411:       MPI_Gather(&mumps->nnz,1,MPIU_INT64,mumps->recvcount,1,MPIU_INT64,0/*master*/,mumps->omp_comm);

1413:       /* Then each computes number of send/recvs */
1414:       if (mumps->is_omp_master) {
1415:         /* Start from 1 since self communication is not done in MPI */
1416:         nreqs = 0;
1417:         for (i=1; i<osize; i++) nreqs += (mumps->recvcount[i]+PETSC_MPI_INT_MAX-1)/PETSC_MPI_INT_MAX;
1418:       } else {
1419:         nreqs = (mumps->nnz+PETSC_MPI_INT_MAX-1)/PETSC_MPI_INT_MAX;
1420:       }
1421:       PetscMalloc1(nreqs*3,&mumps->reqs); /* Triple the requests since we send irn, jcn and val seperately */

1423:       /* The following code is doing a very simple thing: omp_master rank gathers irn/jcn/val from others.
1424:          MPI_Gatherv would be enough if it supports big counts > 2^31-1. Since it does not, and mumps->nnz
1425:          might be a prime number > 2^31-1, we have to slice the message. Note omp_comm_size
1426:          is very small, the current approach should have no extra overhead compared to MPI_Gatherv.
1427:        */
1428:       nreqs = 0; /* counter for actual send/recvs */
1429:       if (mumps->is_omp_master) {
1430:         for (i=0,totnnz=0; i<osize; i++) totnnz += mumps->recvcount[i]; /* totnnz = sum of nnz over omp_comm */
1431:         PetscMalloc2(totnnz,&irn,totnnz,&jcn);
1432:         PetscMalloc1(totnnz,&val);

1434:         /* Self communication */
1435:         PetscArraycpy(irn,mumps->irn,mumps->nnz);
1436:         PetscArraycpy(jcn,mumps->jcn,mumps->nnz);
1437:         PetscArraycpy(val,mumps->val,mumps->nnz);

1439:         /* Replace mumps->irn/jcn etc on master with the newly allocated bigger arrays */
1440:         PetscFree2(mumps->irn,mumps->jcn);
1441:         PetscFree(mumps->val_alloc);
1442:         mumps->nnz = totnnz;
1443:         mumps->irn = irn;
1444:         mumps->jcn = jcn;
1445:         mumps->val = mumps->val_alloc = val;

1447:         irn += mumps->recvcount[0]; /* recvcount[0] is old mumps->nnz on omp rank 0 */
1448:         jcn += mumps->recvcount[0];
1449:         val += mumps->recvcount[0];

1451:         /* Remote communication */
1452:         for (i=1; i<osize; i++) {
1453:           count  = PetscMin(mumps->recvcount[i],PETSC_MPI_INT_MAX);
1454:           remain = mumps->recvcount[i] - count;
1455:           while (count>0) {
1456:             MPI_Irecv(irn,count,MPIU_MUMPSINT,i,mumps->tag,mumps->omp_comm,&mumps->reqs[nreqs++]);
1457:             MPI_Irecv(jcn,count,MPIU_MUMPSINT,i,mumps->tag,mumps->omp_comm,&mumps->reqs[nreqs++]);
1458:             MPI_Irecv(val,count,MPIU_SCALAR,  i,mumps->tag,mumps->omp_comm,&mumps->reqs[nreqs++]);
1459:             irn    += count;
1460:             jcn    += count;
1461:             val    += count;
1462:             count   = PetscMin(remain,PETSC_MPI_INT_MAX);
1463:             remain -= count;
1464:           }
1465:         }
1466:       } else {
1467:         irn    = mumps->irn;
1468:         jcn    = mumps->jcn;
1469:         val    = mumps->val;
1470:         count  = PetscMin(mumps->nnz,PETSC_MPI_INT_MAX);
1471:         remain = mumps->nnz - count;
1472:         while (count>0) {
1473:           MPI_Isend(irn,count,MPIU_MUMPSINT,0,mumps->tag,mumps->omp_comm,&mumps->reqs[nreqs++]);
1474:           MPI_Isend(jcn,count,MPIU_MUMPSINT,0,mumps->tag,mumps->omp_comm,&mumps->reqs[nreqs++]);
1475:           MPI_Isend(val,count,MPIU_SCALAR,  0,mumps->tag,mumps->omp_comm,&mumps->reqs[nreqs++]);
1476:           irn    += count;
1477:           jcn    += count;
1478:           val    += count;
1479:           count   = PetscMin(remain,PETSC_MPI_INT_MAX);
1480:           remain -= count;
1481:         }
1482:       }
1483:     } else {
1484:       nreqs = 0;
1485:       if (mumps->is_omp_master) {
1486:         val = mumps->val + mumps->recvcount[0];
1487:         for (i=1; i<osize; i++) { /* Remote communication only since self data is already in place */
1488:           count  = PetscMin(mumps->recvcount[i],PETSC_MPI_INT_MAX);
1489:           remain = mumps->recvcount[i] - count;
1490:           while (count>0) {
1491:             MPI_Irecv(val,count,MPIU_SCALAR,i,mumps->tag,mumps->omp_comm,&mumps->reqs[nreqs++]);
1492:             val    += count;
1493:             count   = PetscMin(remain,PETSC_MPI_INT_MAX);
1494:             remain -= count;
1495:           }
1496:         }
1497:       } else {
1498:         val    = mumps->val;
1499:         count  = PetscMin(mumps->nnz,PETSC_MPI_INT_MAX);
1500:         remain = mumps->nnz - count;
1501:         while (count>0) {
1502:           MPI_Isend(val,count,MPIU_SCALAR,0,mumps->tag,mumps->omp_comm,&mumps->reqs[nreqs++]);
1503:           val    += count;
1504:           count   = PetscMin(remain,PETSC_MPI_INT_MAX);
1505:           remain -= count;
1506:         }
1507:       }
1508:     }
1509:     MPI_Waitall(nreqs,mumps->reqs,MPI_STATUSES_IGNORE);
1510:     mumps->tag++; /* It is totally fine for above send/recvs to share one mpi tag */
1511:   }
1512:   return(0);
1513: }

1515: PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
1516: {
1517:   Mat_MUMPS      *mumps =(Mat_MUMPS*)(F)->data;
1519:   PetscBool      isMPIAIJ;

1522:   if (mumps->id.INFOG(1) < 0 && !(mumps->id.INFOG(1) == -16 && mumps->id.INFOG(1) == 0)) {
1523:     if (mumps->id.INFOG(1) == -6) {
1524:       PetscInfo2(A,"MatFactorNumeric is called with singular matrix structure, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1525:     }
1526:     PetscInfo2(A,"MatFactorNumeric is called after analysis phase fails, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1527:     return(0);
1528:   }

1530:   (*mumps->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, mumps);
1531:   MatMumpsGatherNonzerosOnMaster(MAT_REUSE_MATRIX,mumps);

1533:   /* numerical factorization phase */
1534:   /*-------------------------------*/
1535:   mumps->id.job = JOB_FACTNUMERIC;
1536:   if (!mumps->id.ICNTL(18)) { /* A is centralized */
1537:     if (!mumps->myid) {
1538:       mumps->id.a = (MumpsScalar*)mumps->val;
1539:     }
1540:   } else {
1541:     mumps->id.a_loc = (MumpsScalar*)mumps->val;
1542:   }
1543:   PetscMUMPS_c(mumps);
1544:   if (mumps->id.INFOG(1) < 0) {
1545:     if (A->erroriffailure) {
1546:       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));
1547:     } else {
1548:       if (mumps->id.INFOG(1) == -10) { /* numerically singular matrix */
1549:         PetscInfo2(F,"matrix is numerically singular, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1550:         F->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1551:       } else if (mumps->id.INFOG(1) == -13) {
1552:         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));
1553:         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1554:       } else if (mumps->id.INFOG(1) == -8 || mumps->id.INFOG(1) == -9 || (-16 < mumps->id.INFOG(1) && mumps->id.INFOG(1) < -10) ) {
1555:         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));
1556:         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1557:       } else {
1558:         PetscInfo2(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1559:         F->factorerrortype = MAT_FACTOR_OTHER;
1560:       }
1561:     }
1562:   }
1563:   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));

1565:   F->assembled    = PETSC_TRUE;
1566:   mumps->matstruc = SAME_NONZERO_PATTERN;
1567:   if (F->schur) { /* reset Schur status to unfactored */
1568: #if defined(PETSC_HAVE_CUDA)
1569:     F->schur->offloadmask = PETSC_OFFLOAD_CPU;
1570: #endif
1571:     if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
1572:       mumps->id.ICNTL(19) = 2;
1573:       MatTranspose(F->schur,MAT_INPLACE_MATRIX,&F->schur);
1574:     }
1575:     MatFactorRestoreSchurComplement(F,NULL,MAT_FACTOR_SCHUR_UNFACTORED);
1576:   }

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

1581:   if (!mumps->is_omp_master) mumps->id.INFO(23) = 0;
1582:   if (mumps->petsc_size > 1) {
1583:     PetscInt    lsol_loc;
1584:     PetscScalar *sol_loc;

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

1588:     /* distributed solution; Create x_seq=sol_loc for repeated use */
1589:     if (mumps->x_seq) {
1590:       VecScatterDestroy(&mumps->scat_sol);
1591:       PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);
1592:       VecDestroy(&mumps->x_seq);
1593:     }
1594:     lsol_loc = mumps->id.INFO(23); /* length of sol_loc */
1595:     PetscMalloc2(lsol_loc,&sol_loc,lsol_loc,&mumps->id.isol_loc);
1596:     mumps->id.lsol_loc = lsol_loc;
1597:     mumps->id.sol_loc = (MumpsScalar*)sol_loc;
1598:     VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq);
1599:   }
1600:   PetscLogFlops(mumps->id.RINFO(2));
1601:   return(0);
1602: }

1604: /* Sets MUMPS options from the options database */
1605: PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A)
1606: {
1607:   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1609:   PetscMUMPSInt  icntl=0;
1610:   PetscInt       info[80],i,ninfo=80;
1611:   PetscBool      flg=PETSC_FALSE;

1614:   PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MUMPS Options","Mat");
1615:   PetscOptionsMUMPSInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);
1616:   if (flg) mumps->id.ICNTL(1) = icntl;
1617:   PetscOptionsMUMPSInt("-mat_mumps_icntl_2","ICNTL(2): output stream for diagnostic printing, statistics, and warning","None",mumps->id.ICNTL(2),&icntl,&flg);
1618:   if (flg) mumps->id.ICNTL(2) = icntl;
1619:   PetscOptionsMUMPSInt("-mat_mumps_icntl_3","ICNTL(3): output stream for global information, collected on the host","None",mumps->id.ICNTL(3),&icntl,&flg);
1620:   if (flg) mumps->id.ICNTL(3) = icntl;

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

1626:   PetscOptionsMUMPSInt("-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);
1627:   if (flg) mumps->id.ICNTL(6) = icntl;

1629:   PetscOptionsMUMPSInt("-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);
1630:   if (flg) {
1631:     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");
1632:     else mumps->id.ICNTL(7) = icntl;
1633:   }

1635:   PetscOptionsMUMPSInt("-mat_mumps_icntl_8","ICNTL(8): scaling strategy (-2 to 8 or 77)","None",mumps->id.ICNTL(8),&mumps->id.ICNTL(8),NULL);
1636:   /* 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() */
1637:   PetscOptionsMUMPSInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",mumps->id.ICNTL(10),&mumps->id.ICNTL(10),NULL);
1638:   PetscOptionsMUMPSInt("-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);
1639:   PetscOptionsMUMPSInt("-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);
1640:   PetscOptionsMUMPSInt("-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);
1641:   PetscOptionsMUMPSInt("-mat_mumps_icntl_14","ICNTL(14): percentage increase in the estimated working space","None",mumps->id.ICNTL(14),&mumps->id.ICNTL(14),NULL);
1642:   PetscOptionsMUMPSInt("-mat_mumps_icntl_19","ICNTL(19): computes the Schur complement","None",mumps->id.ICNTL(19),&mumps->id.ICNTL(19),NULL);
1643:   if (mumps->id.ICNTL(19) <= 0 || mumps->id.ICNTL(19) > 3) { /* reset any schur data (if any) */
1644:     MatDestroy(&F->schur);
1645:     MatMumpsResetSchur_Private(mumps);
1646:   }
1647:   /* PetscOptionsMUMPSInt("-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 */
1648:   /* PetscOptionsMUMPSInt("-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 */

1650:   PetscOptionsMUMPSInt("-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);
1651:   PetscOptionsMUMPSInt("-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);
1652:   PetscOptionsMUMPSInt("-mat_mumps_icntl_24","ICNTL(24): detection of null pivot rows (0 or 1)","None",mumps->id.ICNTL(24),&mumps->id.ICNTL(24),NULL);
1653:   if (mumps->id.ICNTL(24)) {
1654:     mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */
1655:   }

1657:   PetscOptionsMUMPSInt("-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);
1658:   PetscOptionsMUMPSInt("-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);
1659:   PetscOptionsMUMPSInt("-mat_mumps_icntl_27","ICNTL(27): the blocking size for multiple right-hand sides","None",mumps->id.ICNTL(27),&mumps->id.ICNTL(27),NULL);
1660:   PetscOptionsMUMPSInt("-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);
1661:   PetscOptionsMUMPSInt("-mat_mumps_icntl_29","ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis","None",mumps->id.ICNTL(29),&mumps->id.ICNTL(29),NULL);
1662:   /* PetscOptionsMUMPSInt("-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 */
1663:   PetscOptionsMUMPSInt("-mat_mumps_icntl_31","ICNTL(31): indicates which factors may be discarded during factorization","None",mumps->id.ICNTL(31),&mumps->id.ICNTL(31),NULL);
1664:   /* PetscOptionsMUMPSInt("-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 */
1665:   PetscOptionsMUMPSInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),NULL);
1666:   PetscOptionsMUMPSInt("-mat_mumps_icntl_35","ICNTL(35): activates Block Low Rank (BLR) based factorization","None",mumps->id.ICNTL(35),&mumps->id.ICNTL(35),NULL);
1667:   PetscOptionsMUMPSInt("-mat_mumps_icntl_36","ICNTL(36): choice of BLR factorization variant","None",mumps->id.ICNTL(36),&mumps->id.ICNTL(36),NULL);
1668:   PetscOptionsMUMPSInt("-mat_mumps_icntl_38","ICNTL(38): estimated compression rate of LU factors with BLR","None",mumps->id.ICNTL(38),&mumps->id.ICNTL(38),NULL);

1670:   PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),NULL);
1671:   PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",mumps->id.CNTL(2),&mumps->id.CNTL(2),NULL);
1672:   PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),NULL);
1673:   PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",mumps->id.CNTL(4),&mumps->id.CNTL(4),NULL);
1674:   PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",mumps->id.CNTL(5),&mumps->id.CNTL(5),NULL);
1675:   PetscOptionsReal("-mat_mumps_cntl_7","CNTL(7): dropping parameter used during BLR","None",mumps->id.CNTL(7),&mumps->id.CNTL(7),NULL);

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

1679:   PetscOptionsIntArray("-mat_mumps_view_info","request INFO local to each processor","",info,&ninfo,NULL);
1680:   if (ninfo) {
1681:     if (ninfo > 80) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"number of INFO %d must <= 80\n",ninfo);
1682:     PetscMalloc1(ninfo,&mumps->info);
1683:     mumps->ninfo = ninfo;
1684:     for (i=0; i<ninfo; i++) {
1685:       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);
1686:       else  mumps->info[i] = info[i];
1687:     }
1688:   }

1690:   PetscOptionsEnd();
1691:   return(0);
1692: }

1694: PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps)
1695: {
1697:   PetscInt       nthreads=0;

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

1704:   PetscOptionsHasName(NULL,NULL,"-mat_mumps_use_omp_threads",&mumps->use_petsc_omp_support);
1705:   if (mumps->use_petsc_omp_support) nthreads = -1; /* -1 will let PetscOmpCtrlCreate() guess a proper value when user did not supply one */
1706:   PetscOptionsGetInt(NULL,NULL,"-mat_mumps_use_omp_threads",&nthreads,NULL);
1707:   if (mumps->use_petsc_omp_support) {
1708: #if defined(PETSC_HAVE_OPENMP_SUPPORT)
1709:     PetscOmpCtrlCreate(mumps->petsc_comm,nthreads,&mumps->omp_ctrl);
1710:     PetscOmpCtrlGetOmpComms(mumps->omp_ctrl,&mumps->omp_comm,&mumps->mumps_comm,&mumps->is_omp_master);
1711: #else
1712:     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");
1713: #endif
1714:   } else {
1715:     mumps->omp_comm      = PETSC_COMM_SELF;
1716:     mumps->mumps_comm    = mumps->petsc_comm;
1717:     mumps->is_omp_master = PETSC_TRUE;
1718:   }
1719:   MPI_Comm_size(mumps->omp_comm,&mumps->omp_comm_size);
1720:   mumps->reqs = NULL;
1721:   mumps->tag  = 0;

1723:   mumps->id.comm_fortran = MPI_Comm_c2f(mumps->mumps_comm);
1724:   mumps->id.job = JOB_INIT;
1725:   mumps->id.par = 1;  /* host participates factorizaton and solve */
1726:   mumps->id.sym = mumps->sym;

1728:   PetscMUMPS_c(mumps);
1729:   if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in PetscInitializeMUMPS: INFOG(1)=%d\n",mumps->id.INFOG(1));

1731:   /* copy MUMPS default control values from master to slaves. Although slaves do not call MUMPS, they may access these values in code.
1732:      For example, ICNTL(9) is initialized to 1 by MUMPS and slaves check ICNTL(9) in MatSolve_MUMPS.
1733:    */
1734:   MPI_Bcast(mumps->id.icntl,40,MPI_INT,  0,mumps->omp_comm); /* see MUMPS-5.1.2 Manual Section 9 */
1735:   MPI_Bcast(mumps->id.cntl, 15,MPIU_REAL,0,mumps->omp_comm);

1737:   mumps->scat_rhs = NULL;
1738:   mumps->scat_sol = NULL;

1740:   /* set PETSc-MUMPS default options - override MUMPS default */
1741:   mumps->id.ICNTL(3) = 0;
1742:   mumps->id.ICNTL(4) = 0;
1743:   if (mumps->petsc_size == 1) {
1744:     mumps->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
1745:   } else {
1746:     mumps->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
1747:     mumps->id.ICNTL(20) = 0;   /* rhs is in dense format */
1748:     mumps->id.ICNTL(21) = 1;   /* distributed solution */
1749:   }

1751:   /* schur */
1752:   mumps->id.size_schur    = 0;
1753:   mumps->id.listvar_schur = NULL;
1754:   mumps->id.schur         = NULL;
1755:   mumps->sizeredrhs       = 0;
1756:   mumps->schur_sol        = NULL;
1757:   mumps->schur_sizesol    = 0;
1758:   return(0);
1759: }

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

1766:   if (mumps->id.INFOG(1) < 0) {
1767:     if (A->erroriffailure) {
1768:       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
1769:     } else {
1770:       if (mumps->id.INFOG(1) == -6) {
1771:         PetscInfo2(F,"matrix is singular in structure, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1772:         F->factorerrortype = MAT_FACTOR_STRUCT_ZEROPIVOT;
1773:       } else if (mumps->id.INFOG(1) == -5 || mumps->id.INFOG(1) == -7) {
1774:         PetscInfo2(F,"problem of workspace, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1775:         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1776:       } else if (mumps->id.INFOG(1) == -16 && mumps->id.INFOG(1) == 0) {
1777:         PetscInfo(F,"Empty matrix\n");
1778:       } else {
1779:         PetscInfo2(F,"Error reported by MUMPS in analysis phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1780:         F->factorerrortype = MAT_FACTOR_OTHER;
1781:       }
1782:     }
1783:   }
1784:   return(0);
1785: }

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

1796:   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;

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

1801:   (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, mumps);
1802:   MatMumpsGatherNonzerosOnMaster(MAT_INITIAL_MATRIX,mumps);

1804:   /* analysis phase */
1805:   /*----------------*/
1806:   mumps->id.job = JOB_FACTSYMBOLIC;
1807:   mumps->id.n   = M;
1808:   switch (mumps->id.ICNTL(18)) {
1809:   case 0:  /* centralized assembled matrix input */
1810:     if (!mumps->myid) {
1811:       mumps->id.nnz = mumps->nnz;
1812:       mumps->id.irn = mumps->irn;
1813:       mumps->id.jcn = mumps->jcn;
1814:       if (mumps->id.ICNTL(6)>1) mumps->id.a = (MumpsScalar*)mumps->val;
1815:       if (mumps->id.ICNTL(7) == 1) { /* use user-provide matrix ordering - assuming r = c ordering */
1816:         /*
1817:         PetscBool      flag;
1818:         ISEqual(r,c,&flag);
1819:         if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"row_perm != col_perm");
1820:         ISView(r,PETSC_VIEWER_STDOUT_SELF);
1821:          */
1822:         if (!mumps->myid) {
1823:           const PetscInt *idx;
1824:           PetscInt       i;

1826:           PetscMalloc1(M,&mumps->id.perm_in);
1827:           ISGetIndices(r,&idx);
1828:           for (i=0; i<M; i++) {PetscMUMPSIntCast(idx[i]+1,&(mumps->id.perm_in[i]));} /* perm_in[]: start from 1, not 0! */
1829:           ISRestoreIndices(r,&idx);
1830:         }
1831:       }
1832:     }
1833:     break;
1834:   case 3:  /* distributed assembled matrix input (size>1) */
1835:     mumps->id.nnz_loc = mumps->nnz;
1836:     mumps->id.irn_loc = mumps->irn;
1837:     mumps->id.jcn_loc = mumps->jcn;
1838:     if (mumps->id.ICNTL(6)>1) mumps->id.a_loc = (MumpsScalar*)mumps->val;
1839:     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1840:     MatCreateVecs(A,NULL,&b);
1841:     VecScatterCreateToZero(b,&mumps->scat_rhs,&mumps->b_seq);
1842:     VecDestroy(&b);
1843:     break;
1844:   }
1845:   PetscMUMPS_c(mumps);
1846:   MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);

1848:   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1849:   F->ops->solve           = MatSolve_MUMPS;
1850:   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1851:   F->ops->matsolve        = MatMatSolve_MUMPS;
1852:   F->ops->mattransposesolve = MatMatTransposeSolve_MUMPS;
1853:   return(0);
1854: }

1856: /* Note the Petsc r and c permutations are ignored */
1857: PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1858: {
1859:   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1861:   Vec            b;
1862:   const PetscInt M = A->rmap->N;

1865:   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;

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

1870:   (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, mumps);
1871:   MatMumpsGatherNonzerosOnMaster(MAT_INITIAL_MATRIX,mumps);

1873:   /* analysis phase */
1874:   /*----------------*/
1875:   mumps->id.job = JOB_FACTSYMBOLIC;
1876:   mumps->id.n   = M;
1877:   switch (mumps->id.ICNTL(18)) {
1878:   case 0:  /* centralized assembled matrix input */
1879:     if (!mumps->myid) {
1880:       mumps->id.nnz = mumps->nnz;
1881:       mumps->id.irn = mumps->irn;
1882:       mumps->id.jcn = mumps->jcn;
1883:       if (mumps->id.ICNTL(6)>1) {
1884:         mumps->id.a = (MumpsScalar*)mumps->val;
1885:       }
1886:     }
1887:     break;
1888:   case 3:  /* distributed assembled matrix input (size>1) */
1889:     mumps->id.nnz_loc = mumps->nnz;
1890:     mumps->id.irn_loc = mumps->irn;
1891:     mumps->id.jcn_loc = mumps->jcn;
1892:     if (mumps->id.ICNTL(6)>1) {
1893:       mumps->id.a_loc = (MumpsScalar*)mumps->val;
1894:     }
1895:     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1896:     MatCreateVecs(A,NULL,&b);
1897:     VecScatterCreateToZero(b,&mumps->scat_rhs,&mumps->b_seq);
1898:     VecDestroy(&b);
1899:     break;
1900:   }
1901:   PetscMUMPS_c(mumps);
1902:   MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);

1904:   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1905:   F->ops->solve           = MatSolve_MUMPS;
1906:   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1907:   return(0);
1908: }

1910: /* Note the Petsc r permutation and factor info are ignored */
1911: PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
1912: {
1913:   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1915:   Vec            b;
1916:   const PetscInt M = A->rmap->N;

1919:   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;

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

1924:   (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, mumps);
1925:   MatMumpsGatherNonzerosOnMaster(MAT_INITIAL_MATRIX,mumps);

1927:   /* analysis phase */
1928:   /*----------------*/
1929:   mumps->id.job = JOB_FACTSYMBOLIC;
1930:   mumps->id.n   = M;
1931:   switch (mumps->id.ICNTL(18)) {
1932:   case 0:  /* centralized assembled matrix input */
1933:     if (!mumps->myid) {
1934:       mumps->id.nnz = mumps->nnz;
1935:       mumps->id.irn = mumps->irn;
1936:       mumps->id.jcn = mumps->jcn;
1937:       if (mumps->id.ICNTL(6)>1) {
1938:         mumps->id.a = (MumpsScalar*)mumps->val;
1939:       }
1940:     }
1941:     break;
1942:   case 3:  /* distributed assembled matrix input (size>1) */
1943:     mumps->id.nnz_loc = mumps->nnz;
1944:     mumps->id.irn_loc = mumps->irn;
1945:     mumps->id.jcn_loc = mumps->jcn;
1946:     if (mumps->id.ICNTL(6)>1) {
1947:       mumps->id.a_loc = (MumpsScalar*)mumps->val;
1948:     }
1949:     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1950:     MatCreateVecs(A,NULL,&b);
1951:     VecScatterCreateToZero(b,&mumps->scat_rhs,&mumps->b_seq);
1952:     VecDestroy(&b);
1953:     break;
1954:   }
1955:   PetscMUMPS_c(mumps);
1956:   MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);

1958:   F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
1959:   F->ops->solve                 = MatSolve_MUMPS;
1960:   F->ops->solvetranspose        = MatSolve_MUMPS;
1961:   F->ops->matsolve              = MatMatSolve_MUMPS;
1962:   F->ops->mattransposesolve     = MatMatTransposeSolve_MUMPS;
1963: #if defined(PETSC_USE_COMPLEX)
1964:   F->ops->getinertia = NULL;
1965: #else
1966:   F->ops->getinertia = MatGetInertia_SBAIJMUMPS;
1967: #endif
1968:   return(0);
1969: }

1971: PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
1972: {
1973:   PetscErrorCode    ierr;
1974:   PetscBool         iascii;
1975:   PetscViewerFormat format;
1976:   Mat_MUMPS         *mumps=(Mat_MUMPS*)A->data;

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

1982:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1983:   if (iascii) {
1984:     PetscViewerGetFormat(viewer,&format);
1985:     if (format == PETSC_VIEWER_ASCII_INFO) {
1986:       PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");
1987:       PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                   %d \n",mumps->id.sym);
1988:       PetscViewerASCIIPrintf(viewer,"  PAR (host participation):            %d \n",mumps->id.par);
1989:       PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):         %d \n",mumps->id.ICNTL(1));
1990:       PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg): %d \n",mumps->id.ICNTL(2));
1991:       PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):   %d \n",mumps->id.ICNTL(3));
1992:       PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):        %d \n",mumps->id.ICNTL(4));
1993:       PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):         %d \n",mumps->id.ICNTL(5));
1994:       PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):        %d \n",mumps->id.ICNTL(6));
1995:       PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (sequential matrix ordering):%d \n",mumps->id.ICNTL(7));
1996:       PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scaling strategy):        %d \n",mumps->id.ICNTL(8));
1997:       PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements):  %d \n",mumps->id.ICNTL(10));
1998:       PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):          %d \n",mumps->id.ICNTL(11));
1999:       if (mumps->id.ICNTL(11)>0) {
2000:         PetscViewerASCIIPrintf(viewer,"    RINFOG(4) (inf norm of input mat):        %g\n",mumps->id.RINFOG(4));
2001:         PetscViewerASCIIPrintf(viewer,"    RINFOG(5) (inf norm of solution):         %g\n",mumps->id.RINFOG(5));
2002:         PetscViewerASCIIPrintf(viewer,"    RINFOG(6) (inf norm of residual):         %g\n",mumps->id.RINFOG(6));
2003:         PetscViewerASCIIPrintf(viewer,"    RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8));
2004:         PetscViewerASCIIPrintf(viewer,"    RINFOG(9) (error estimate):               %g \n",mumps->id.RINFOG(9));
2005:         PetscViewerASCIIPrintf(viewer,"    RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11));
2006:       }
2007:       PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",mumps->id.ICNTL(12));
2008:       PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (sequential factorization of the root node):  %d \n",mumps->id.ICNTL(13));
2009:       PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",mumps->id.ICNTL(14));
2010:       /* ICNTL(15-17) not used */
2011:       PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",mumps->id.ICNTL(18));
2012:       PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Schur complement info):                      %d \n",mumps->id.ICNTL(19));
2013:       PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",mumps->id.ICNTL(20));
2014:       PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d \n",mumps->id.ICNTL(21));
2015:       PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",mumps->id.ICNTL(22));
2016:       PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",mumps->id.ICNTL(23));

2018:       PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",mumps->id.ICNTL(24));
2019:       PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",mumps->id.ICNTL(25));
2020:       PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",mumps->id.ICNTL(26));
2021:       PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",mumps->id.ICNTL(27));
2022:       PetscViewerASCIIPrintf(viewer,"  ICNTL(28) (use parallel or sequential ordering):        %d \n",mumps->id.ICNTL(28));
2023:       PetscViewerASCIIPrintf(viewer,"  ICNTL(29) (parallel ordering):                          %d \n",mumps->id.ICNTL(29));

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

2032:       PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",mumps->id.CNTL(1));
2033:       PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2));
2034:       PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",mumps->id.CNTL(3));
2035:       PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",mumps->id.CNTL(4));
2036:       PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",mumps->id.CNTL(5));
2037:       PetscViewerASCIIPrintf(viewer,"  CNTL(7) (dropping parameter for BLR):       %g \n",mumps->id.CNTL(7));

2039:       /* infomation local to each processor */
2040:       PetscViewerASCIIPrintf(viewer, "  RINFO(1) (local estimated flops for the elimination after analysis): \n");
2041:       PetscViewerASCIIPushSynchronized(viewer);
2042:       PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %g \n",mumps->myid,mumps->id.RINFO(1));
2043:       PetscViewerFlush(viewer);
2044:       PetscViewerASCIIPrintf(viewer, "  RINFO(2) (local estimated flops for the assembly after factorization): \n");
2045:       PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(2));
2046:       PetscViewerFlush(viewer);
2047:       PetscViewerASCIIPrintf(viewer, "  RINFO(3) (local estimated flops for the elimination after factorization): \n");
2048:       PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(3));
2049:       PetscViewerFlush(viewer);

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

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

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

2063:       if (mumps->ninfo && mumps->ninfo <= 80){
2064:         PetscInt i;
2065:         for (i=0; i<mumps->ninfo; i++){
2066:           PetscViewerASCIIPrintf(viewer, "  INFO(%d): \n",mumps->info[i]);
2067:           PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(mumps->info[i]));
2068:           PetscViewerFlush(viewer);
2069:         }
2070:       }
2071:       PetscViewerASCIIPopSynchronized(viewer);

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

2079:         PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(3));
2080:         PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(4));
2081:         PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",mumps->id.INFOG(5));
2082:         PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",mumps->id.INFOG(6));
2083:         PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively use after analysis): %d \n",mumps->id.INFOG(7));
2084:         PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",mumps->id.INFOG(8));
2085:         PetscViewerASCIIPrintf(viewer,"  INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",mumps->id.INFOG(9));
2086:         PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",mumps->id.INFOG(10));
2087:         PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",mumps->id.INFOG(11));
2088:         PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",mumps->id.INFOG(12));
2089:         PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",mumps->id.INFOG(13));
2090:         PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",mumps->id.INFOG(14));
2091:         PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",mumps->id.INFOG(15));
2092:         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));
2093:         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));
2094:         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));
2095:         PetscViewerASCIIPrintf(viewer,"  INFOG(19) (size of all MUMPS internal data allocated during factorization: sum over all processors): %d \n",mumps->id.INFOG(19));
2096:         PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",mumps->id.INFOG(20));
2097:         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));
2098:         PetscViewerASCIIPrintf(viewer,"  INFOG(22) (size in MB of memory effectively used during factorization - sum over all processors): %d \n",mumps->id.INFOG(22));
2099:         PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",mumps->id.INFOG(23));
2100:         PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",mumps->id.INFOG(24));
2101:         PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",mumps->id.INFOG(25));
2102:         PetscViewerASCIIPrintf(viewer,"  INFOG(28) (after factorization: number of null pivots encountered): %d\n",mumps->id.INFOG(28));
2103:         PetscViewerASCIIPrintf(viewer,"  INFOG(29) (after factorization: effective number of entries in the factors (sum over all processors)): %d\n",mumps->id.INFOG(29));
2104:         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));
2105:         PetscViewerASCIIPrintf(viewer,"  INFOG(32) (after analysis: type of analysis done): %d\n",mumps->id.INFOG(32));
2106:         PetscViewerASCIIPrintf(viewer,"  INFOG(33) (value used for ICNTL(8)): %d\n",mumps->id.INFOG(33));
2107:         PetscViewerASCIIPrintf(viewer,"  INFOG(34) (exponent of the determinant if determinant is requested): %d\n",mumps->id.INFOG(34));
2108:         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));
2109:         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));
2110:         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));
2111:         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));
2112:         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));
2113:       }
2114:     }
2115:   }
2116:   return(0);
2117: }

2119: PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
2120: {
2121:   Mat_MUMPS *mumps =(Mat_MUMPS*)A->data;

2124:   info->block_size        = 1.0;
2125:   info->nz_allocated      = mumps->id.INFOG(20);
2126:   info->nz_used           = mumps->id.INFOG(20);
2127:   info->nz_unneeded       = 0.0;
2128:   info->assemblies        = 0.0;
2129:   info->mallocs           = 0.0;
2130:   info->memory            = 0.0;
2131:   info->fill_ratio_given  = 0;
2132:   info->fill_ratio_needed = 0;
2133:   info->factor_mallocs    = 0;
2134:   return(0);
2135: }

2137: /* -------------------------------------------------------------------------------------------*/
2138: PetscErrorCode MatFactorSetSchurIS_MUMPS(Mat F, IS is)
2139: {
2140:   Mat_MUMPS         *mumps =(Mat_MUMPS*)F->data;
2141:   const PetscScalar *arr;
2142:   const PetscInt    *idxs;
2143:   PetscInt          size,i;
2144:   PetscErrorCode    ierr;

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

2151:     ls   = mumps->myid ? (size ? PETSC_FALSE : PETSC_TRUE) : PETSC_TRUE; /* always true on root; false on others if their size != 0 */
2152:     MPI_Allreduce(&ls,&gs,1,MPIU_BOOL,MPI_LAND,mumps->petsc_comm);
2153:     if (!gs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MUMPS distributed parallel Schur complements not yet supported from PETSc\n");
2154:   }

2156:   /* Schur complement matrix */
2157:   MatDestroy(&F->schur);
2158:   MatCreateSeqDense(PETSC_COMM_SELF,size,size,NULL,&F->schur);
2159:   MatDenseGetArrayRead(F->schur,&arr);
2160:   mumps->id.schur      = (MumpsScalar*)arr;
2161:   mumps->id.size_schur = size;
2162:   mumps->id.schur_lld  = size;
2163:   MatDenseRestoreArrayRead(F->schur,&arr);
2164:   if (mumps->sym == 1) {
2165:     MatSetOption(F->schur,MAT_SPD,PETSC_TRUE);
2166:   }

2168:   /* MUMPS expects Fortran style indices */
2169:   PetscFree(mumps->id.listvar_schur);
2170:   PetscMalloc1(size,&mumps->id.listvar_schur);
2171:   ISGetIndices(is,&idxs);
2172:   for (i=0; i<size; i++) {PetscMUMPSIntCast(idxs[i]+1,&(mumps->id.listvar_schur[i]));}
2173:   ISRestoreIndices(is,&idxs);
2174:   if (mumps->petsc_size > 1) {
2175:     mumps->id.ICNTL(19) = 1; /* MUMPS returns Schur centralized on the host */
2176:   } else {
2177:     if (F->factortype == MAT_FACTOR_LU) {
2178:       mumps->id.ICNTL(19) = 3; /* MUMPS returns full matrix */
2179:     } else {
2180:       mumps->id.ICNTL(19) = 2; /* MUMPS returns lower triangular part */
2181:     }
2182:   }
2183:   /* set a special value of ICNTL (not handled my MUMPS) to be used in the solve phase by PETSc */
2184:   mumps->id.ICNTL(26) = -1;
2185:   return(0);
2186: }

2188: /* -------------------------------------------------------------------------------------------*/
2189: PetscErrorCode MatFactorCreateSchurComplement_MUMPS(Mat F,Mat* S)
2190: {
2191:   Mat            St;
2192:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
2193:   PetscScalar    *array;
2194: #if defined(PETSC_USE_COMPLEX)
2195:   PetscScalar    im = PetscSqrtScalar((PetscScalar)-1.0);
2196: #endif

2200:   if (!mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it");
2201:   MatCreate(PETSC_COMM_SELF,&St);
2202:   MatSetSizes(St,PETSC_DECIDE,PETSC_DECIDE,mumps->id.size_schur,mumps->id.size_schur);
2203:   MatSetType(St,MATDENSE);
2204:   MatSetUp(St);
2205:   MatDenseGetArray(St,&array);
2206:   if (!mumps->sym) { /* MUMPS always return a full matrix */
2207:     if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
2208:       PetscInt i,j,N=mumps->id.size_schur;
2209:       for (i=0;i<N;i++) {
2210:         for (j=0;j<N;j++) {
2211: #if !defined(PETSC_USE_COMPLEX)
2212:           PetscScalar val = mumps->id.schur[i*N+j];
2213: #else
2214:           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
2215: #endif
2216:           array[j*N+i] = val;
2217:         }
2218:       }
2219:     } else { /* stored by columns */
2220:       PetscArraycpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur);
2221:     }
2222:   } else { /* either full or lower-triangular (not packed) */
2223:     if (mumps->id.ICNTL(19) == 2) { /* lower triangular stored by columns */
2224:       PetscInt i,j,N=mumps->id.size_schur;
2225:       for (i=0;i<N;i++) {
2226:         for (j=i;j<N;j++) {
2227: #if !defined(PETSC_USE_COMPLEX)
2228:           PetscScalar val = mumps->id.schur[i*N+j];
2229: #else
2230:           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
2231: #endif
2232:           array[i*N+j] = val;
2233:           array[j*N+i] = val;
2234:         }
2235:       }
2236:     } else if (mumps->id.ICNTL(19) == 3) { /* full matrix */
2237:       PetscArraycpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur);
2238:     } else { /* ICNTL(19) == 1 lower triangular stored by rows */
2239:       PetscInt i,j,N=mumps->id.size_schur;
2240:       for (i=0;i<N;i++) {
2241:         for (j=0;j<i+1;j++) {
2242: #if !defined(PETSC_USE_COMPLEX)
2243:           PetscScalar val = mumps->id.schur[i*N+j];
2244: #else
2245:           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
2246: #endif
2247:           array[i*N+j] = val;
2248:           array[j*N+i] = val;
2249:         }
2250:       }
2251:     }
2252:   }
2253:   MatDenseRestoreArray(St,&array);
2254:   *S   = St;
2255:   return(0);
2256: }

2258: /* -------------------------------------------------------------------------------------------*/
2259: PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival)
2260: {
2262:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

2265:   PetscMUMPSIntCast(ival,&mumps->id.ICNTL(icntl));
2266:   return(0);
2267: }

2269: PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt *ival)
2270: {
2271:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

2274:   *ival = mumps->id.ICNTL(icntl);
2275:   return(0);
2276: }

2278: /*@
2279:   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()

2281:    Logically Collective on Mat

2283:    Input Parameters:
2284: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2285: .  icntl - index of MUMPS parameter array ICNTL()
2286: -  ival - value of MUMPS ICNTL(icntl)

2288:   Options Database:
2289: .   -mat_mumps_icntl_<icntl> <ival>

2291:    Level: beginner

2293:    References:
2294: .     MUMPS Users' Guide

2296: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2297:  @*/
2298: PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival)
2299: {

2304:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2307:   PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));
2308:   return(0);
2309: }

2311: /*@
2312:   MatMumpsGetIcntl - Get MUMPS parameter ICNTL()

2314:    Logically Collective on Mat

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

2320:   Output Parameter:
2321: .  ival - value of MUMPS ICNTL(icntl)

2323:    Level: beginner

2325:    References:
2326: .     MUMPS Users' Guide

2328: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2329: @*/
2330: PetscErrorCode MatMumpsGetIcntl(Mat F,PetscInt icntl,PetscInt *ival)
2331: {

2336:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2339:   PetscUseMethod(F,"MatMumpsGetIcntl_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));
2340:   return(0);
2341: }

2343: /* -------------------------------------------------------------------------------------------*/
2344: PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal val)
2345: {
2346:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

2349:   mumps->id.CNTL(icntl) = val;
2350:   return(0);
2351: }

2353: PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal *val)
2354: {
2355:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

2358:   *val = mumps->id.CNTL(icntl);
2359:   return(0);
2360: }

2362: /*@
2363:   MatMumpsSetCntl - Set MUMPS parameter CNTL()

2365:    Logically Collective on Mat

2367:    Input Parameters:
2368: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2369: .  icntl - index of MUMPS parameter array CNTL()
2370: -  val - value of MUMPS CNTL(icntl)

2372:   Options Database:
2373: .   -mat_mumps_cntl_<icntl> <val>

2375:    Level: beginner

2377:    References:
2378: .     MUMPS Users' Guide

2380: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2381: @*/
2382: PetscErrorCode MatMumpsSetCntl(Mat F,PetscInt icntl,PetscReal val)
2383: {

2388:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2391:   PetscTryMethod(F,"MatMumpsSetCntl_C",(Mat,PetscInt,PetscReal),(F,icntl,val));
2392:   return(0);
2393: }

2395: /*@
2396:   MatMumpsGetCntl - Get MUMPS parameter CNTL()

2398:    Logically Collective on Mat

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

2404:   Output Parameter:
2405: .  val - value of MUMPS CNTL(icntl)

2407:    Level: beginner

2409:    References:
2410: .      MUMPS Users' Guide

2412: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2413: @*/
2414: PetscErrorCode MatMumpsGetCntl(Mat F,PetscInt icntl,PetscReal *val)
2415: {

2420:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2423:   PetscUseMethod(F,"MatMumpsGetCntl_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));
2424:   return(0);
2425: }

2427: PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F,PetscInt icntl,PetscInt *info)
2428: {
2429:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

2432:   *info = mumps->id.INFO(icntl);
2433:   return(0);
2434: }

2436: PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F,PetscInt icntl,PetscInt *infog)
2437: {
2438:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

2441:   *infog = mumps->id.INFOG(icntl);
2442:   return(0);
2443: }

2445: PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfo)
2446: {
2447:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

2450:   *rinfo = mumps->id.RINFO(icntl);
2451:   return(0);
2452: }

2454: PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfog)
2455: {
2456:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

2459:   *rinfog = mumps->id.RINFOG(icntl);
2460:   return(0);
2461: }

2463: PetscErrorCode MatMumpsGetInverse_MUMPS(Mat F,Mat spRHS)
2464: {
2466:   Mat            Bt = NULL,Btseq = NULL;
2467:   PetscBool      flg;
2468:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
2469:   PetscScalar    *aa;
2470:   PetscInt       spnr,*ia,*ja,M,nrhs;

2474:   PetscObjectTypeCompare((PetscObject)spRHS,MATTRANSPOSEMAT,&flg);
2475:   if (flg) {
2476:     MatTransposeGetMat(spRHS,&Bt);
2477:   } else SETERRQ(PetscObjectComm((PetscObject)spRHS),PETSC_ERR_ARG_WRONG,"Matrix spRHS must be type MATTRANSPOSEMAT matrix");

2479:   MatMumpsSetIcntl(F,30,1);

2481:   if (mumps->petsc_size > 1) {
2482:     Mat_MPIAIJ *b = (Mat_MPIAIJ*)Bt->data;
2483:     Btseq = b->A;
2484:   } else {
2485:     Btseq = Bt;
2486:   }

2488:   MatGetSize(spRHS,&M,&nrhs);
2489:   mumps->id.nrhs = nrhs;
2490:   mumps->id.lrhs = M;
2491:   mumps->id.rhs  = NULL;

2493:   if (!mumps->myid) {
2494:     MatSeqAIJGetArray(Btseq,&aa);
2495:     MatGetRowIJ(Btseq,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&flg);
2496:     if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot get IJ structure");
2497:     PetscMUMPSIntCSRCast(mumps,spnr,ia,ja,&mumps->id.irhs_ptr,&mumps->id.irhs_sparse,&mumps->id.nz_rhs);
2498:     mumps->id.rhs_sparse  = (MumpsScalar*)aa;
2499:   } else {
2500:     mumps->id.irhs_ptr    = NULL;
2501:     mumps->id.irhs_sparse = NULL;
2502:     mumps->id.nz_rhs      = 0;
2503:     mumps->id.rhs_sparse  = NULL;
2504:   }
2505:   mumps->id.ICNTL(20)   = 1; /* rhs is sparse */
2506:   mumps->id.ICNTL(21)   = 0; /* solution is in assembled centralized format */

2508:   /* solve phase */
2509:   /*-------------*/
2510:   mumps->id.job = JOB_SOLVE;
2511:   PetscMUMPS_c(mumps);
2512:   if (mumps->id.INFOG(1) < 0)
2513:     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));

2515:   if (!mumps->myid) {
2516:     MatSeqAIJRestoreArray(Btseq,&aa);
2517:     MatRestoreRowIJ(Btseq,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&flg);
2518:     if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot get IJ structure");
2519:   }
2520:   return(0);
2521: }

2523: /*@
2524:   MatMumpsGetInverse - Get user-specified set of entries in inverse of A

2526:    Logically Collective on Mat

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

2532:   Output Parameter:
2533: . spRHS - requested entries of inverse of A

2535:    Level: beginner

2537:    References:
2538: .      MUMPS Users' Guide

2540: .seealso: MatGetFactor(), MatCreateTranspose()
2541: @*/
2542: PetscErrorCode MatMumpsGetInverse(Mat F,Mat spRHS)
2543: {

2548:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2549:   PetscUseMethod(F,"MatMumpsGetInverse_C",(Mat,Mat),(F,spRHS));
2550:   return(0);
2551: }

2553: PetscErrorCode MatMumpsGetInverseTranspose_MUMPS(Mat F,Mat spRHST)
2554: {
2556:   Mat            spRHS;

2559:   MatCreateTranspose(spRHST,&spRHS);
2560:   MatMumpsGetInverse_MUMPS(F,spRHS);
2561:   MatDestroy(&spRHS);
2562:   return(0);
2563: }

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

2568:    Logically Collective on Mat

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

2574:   Output Parameter:
2575: . spRHST - requested entries of inverse of A^T

2577:    Level: beginner

2579:    References:
2580: .      MUMPS Users' Guide

2582: .seealso: MatGetFactor(), MatCreateTranspose(), MatMumpsGetInverse()
2583: @*/
2584: PetscErrorCode MatMumpsGetInverseTranspose(Mat F,Mat spRHST)
2585: {
2587:   PetscBool      flg;

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

2595:   PetscUseMethod(F,"MatMumpsGetInverseTranspose_C",(Mat,Mat),(F,spRHST));
2596:   return(0);
2597: }

2599: /*@
2600:   MatMumpsGetInfo - Get MUMPS parameter INFO()

2602:    Logically Collective on Mat

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

2608:   Output Parameter:
2609: .  ival - value of MUMPS INFO(icntl)

2611:    Level: beginner

2613:    References:
2614: .      MUMPS Users' Guide

2616: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2617: @*/
2618: PetscErrorCode MatMumpsGetInfo(Mat F,PetscInt icntl,PetscInt *ival)
2619: {

2624:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2626:   PetscUseMethod(F,"MatMumpsGetInfo_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));
2627:   return(0);
2628: }

2630: /*@
2631:   MatMumpsGetInfog - Get MUMPS parameter INFOG()

2633:    Logically Collective on Mat

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

2639:   Output Parameter:
2640: .  ival - value of MUMPS INFOG(icntl)

2642:    Level: beginner

2644:    References:
2645: .      MUMPS Users' Guide

2647: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2648: @*/
2649: PetscErrorCode MatMumpsGetInfog(Mat F,PetscInt icntl,PetscInt *ival)
2650: {

2655:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2657:   PetscUseMethod(F,"MatMumpsGetInfog_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));
2658:   return(0);
2659: }

2661: /*@
2662:   MatMumpsGetRinfo - Get MUMPS parameter RINFO()

2664:    Logically Collective on Mat

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

2670:   Output Parameter:
2671: .  val - value of MUMPS RINFO(icntl)

2673:    Level: beginner

2675:    References:
2676: .       MUMPS Users' Guide

2678: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2679: @*/
2680: PetscErrorCode MatMumpsGetRinfo(Mat F,PetscInt icntl,PetscReal *val)
2681: {

2686:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2688:   PetscUseMethod(F,"MatMumpsGetRinfo_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));
2689:   return(0);
2690: }

2692: /*@
2693:   MatMumpsGetRinfog - Get MUMPS parameter RINFOG()

2695:    Logically Collective on Mat

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

2701:   Output Parameter:
2702: .  val - value of MUMPS RINFOG(icntl)

2704:    Level: beginner

2706:    References:
2707: .      MUMPS Users' Guide

2709: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2710: @*/
2711: PetscErrorCode MatMumpsGetRinfog(Mat F,PetscInt icntl,PetscReal *val)
2712: {

2717:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2719:   PetscUseMethod(F,"MatMumpsGetRinfog_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));
2720:   return(0);
2721: }

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

2727:   Works with MATAIJ and MATSBAIJ matrices

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

2731:   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.

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

2735:   Options Database Keys:
2736: +  -mat_mumps_icntl_1 - ICNTL(1): output stream for error messages
2737: .  -mat_mumps_icntl_2 - ICNTL(2): output stream for diagnostic printing, statistics, and warning
2738: .  -mat_mumps_icntl_3 -  ICNTL(3): output stream for global information, collected on the host
2739: .  -mat_mumps_icntl_4 -  ICNTL(4): level of printing (0 to 4)
2740: .  -mat_mumps_icntl_6 - ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7)
2741: .  -mat_mumps_icntl_7 - ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis
2742: .  -mat_mumps_icntl_8  - ICNTL(8): scaling strategy (-2 to 8 or 77)
2743: .  -mat_mumps_icntl_10  - ICNTL(10): max num of refinements
2744: .  -mat_mumps_icntl_11  - ICNTL(11): statistics related to an error analysis (via -ksp_view)
2745: .  -mat_mumps_icntl_12  - ICNTL(12): an ordering strategy for symmetric matrices (0 to 3)
2746: .  -mat_mumps_icntl_13  - ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting
2747: .  -mat_mumps_icntl_14  - ICNTL(14): percentage increase in the estimated working space
2748: .  -mat_mumps_icntl_19  - ICNTL(19): computes the Schur complement
2749: .  -mat_mumps_icntl_22  - ICNTL(22): in-core/out-of-core factorization and solve (0 or 1)
2750: .  -mat_mumps_icntl_23  - ICNTL(23): max size of the working memory (MB) that can allocate per processor
2751: .  -mat_mumps_icntl_24  - ICNTL(24): detection of null pivot rows (0 or 1)
2752: .  -mat_mumps_icntl_25  - ICNTL(25): compute a solution of a deficient matrix and a null space basis
2753: .  -mat_mumps_icntl_26  - ICNTL(26): drives the solution phase if a Schur complement matrix
2754: .  -mat_mumps_icntl_28  - ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering
2755: .  -mat_mumps_icntl_29 - ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis
2756: .  -mat_mumps_icntl_30 - ICNTL(30): compute user-specified set of entries in inv(A)
2757: .  -mat_mumps_icntl_31 - ICNTL(31): indicates which factors may be discarded during factorization
2758: .  -mat_mumps_icntl_33 - ICNTL(33): compute determinant
2759: .  -mat_mumps_icntl_35 - ICNTL(35): level of activation of BLR (Block Low-Rank) feature
2760: .  -mat_mumps_icntl_36 - ICNTL(36): controls the choice of BLR factorization variant
2761: .  -mat_mumps_icntl_38 - ICNTL(38): sets the estimated compression rate of LU factors with BLR
2762: .  -mat_mumps_cntl_1  - CNTL(1): relative pivoting threshold
2763: .  -mat_mumps_cntl_2  -  CNTL(2): stopping criterion of refinement
2764: .  -mat_mumps_cntl_3 - CNTL(3): absolute pivoting threshold
2765: .  -mat_mumps_cntl_4 - CNTL(4): value for static pivoting
2766: .  -mat_mumps_cntl_5 - CNTL(5): fixation for null pivots
2767: .  -mat_mumps_cntl_7 - CNTL(7): precision of the dropping parameter used during BLR factorization
2768: -  -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.
2769:                                    Default might be the number of cores per CPU package (socket) as reported by hwloc and suggested by the MUMPS manual.

2771:   Level: beginner

2773:     Notes:
2774:     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.

2776:     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
2777: $          KSPGetPC(ksp,&pc);
2778: $          PCFactorGetMatrix(pc,&mat);
2779: $          MatMumpsGetInfo(mat,....);
2780: $          MatMumpsGetInfog(mat,....); etc.
2781:            Or you can run with -ksp_error_if_not_converged and the program will be stopped and the information printed in the error message.

2783:    Two modes to run MUMPS/PETSc with OpenMP

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

2788: $     -mat_mumps_use_omp_threads [m] and run your code with as many MPI ranks as the number of cores. For example,
2789: $     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"

2791:    To run MUMPS in MPI+OpenMP hybrid mode (i.e., enable multithreading in MUMPS), but still run the non-MUMPS part
2792:    (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
2793:    (or --with-hwloc), and have an MPI that supports MPI-3.0's process shared memory (which is usually available). Since MUMPS calls BLAS
2794:    libraries, to really get performance, you should have multithreaded BLAS libraries such as Intel MKL, AMD ACML, Cray libSci or OpenBLAS
2795:    (PETSc will automatically try to utilized a threaded BLAS if --with-openmp is provided).

2797:    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
2798:    processes on each compute node. Listing the processes in rank ascending order, we split processes on a node into consecutive groups of
2799:    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
2800:    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
2801:    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.
2802:    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,
2803:    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
2804:    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
2805:    cores in socket 1, that definitely hurts locality. On the other hand, if you map MPI ranks consecutively on the two sockets, then the
2806:    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.
2807:    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
2808:    examine the mapping result.

2810:    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,
2811:    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
2812:    calls omp_set_num_threads(m) internally before calling MUMPS.

2814:    References:
2815: +   1. - Heroux, Michael A., R. Brightwell, and Michael M. Wolf. "Bi-modal MPI and MPI+ threads computing on scalable multicore systems." IJHPCA (Submitted) (2011).
2816: -   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.

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

2820: M*/

2822: static PetscErrorCode MatFactorGetSolverType_mumps(Mat A,MatSolverType *type)
2823: {
2825:   *type = MATSOLVERMUMPS;
2826:   return(0);
2827: }

2829: /* MatGetFactor for Seq and MPI AIJ matrices */
2830: static PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
2831: {
2832:   Mat            B;
2834:   Mat_MUMPS      *mumps;
2835:   PetscBool      isSeqAIJ;

2838:  #if defined(PETSC_USE_COMPLEX)
2839:   if (A->hermitian && !A->symmetric && ftype == MAT_FACTOR_CHOLESKY) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Hermitian CHOLESKY Factor is not supported");
2840:  #endif
2841:   /* Create the factorization matrix */
2842:   PetscObjectBaseTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);
2843:   MatCreate(PetscObjectComm((PetscObject)A),&B);
2844:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
2845:   PetscStrallocpy("mumps",&((PetscObject)B)->type_name);
2846:   MatSetUp(B);

2848:   PetscNewLog(B,&mumps);

2850:   B->ops->view    = MatView_MUMPS;
2851:   B->ops->getinfo = MatGetInfo_MUMPS;

2853:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);
2854:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);
2855:   PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);
2856:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);
2857:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);
2858:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);
2859:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);
2860:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);
2861:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);
2862:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);
2863:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);
2864:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverse_C",MatMumpsGetInverse_MUMPS);
2865:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverseTranspose_C",MatMumpsGetInverseTranspose_MUMPS);

2867:   if (ftype == MAT_FACTOR_LU) {
2868:     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
2869:     B->factortype            = MAT_FACTOR_LU;
2870:     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
2871:     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
2872:     mumps->sym = 0;
2873:   } else {
2874:     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2875:     B->factortype                  = MAT_FACTOR_CHOLESKY;
2876:     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
2877:     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
2878: #if defined(PETSC_USE_COMPLEX)
2879:     mumps->sym = 2;
2880: #else
2881:     if (A->spd_set && A->spd) mumps->sym = 1;
2882:     else                      mumps->sym = 2;
2883: #endif
2884:   }

2886:   /* set solvertype */
2887:   PetscFree(B->solvertype);
2888:   PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);

2890:   B->ops->destroy = MatDestroy_MUMPS;
2891:   B->data         = (void*)mumps;

2893:   PetscInitializeMUMPS(A,mumps);

2895:   *F = B;
2896:   return(0);
2897: }

2899: /* MatGetFactor for Seq and MPI SBAIJ matrices */
2900: static PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
2901: {
2902:   Mat            B;
2904:   Mat_MUMPS      *mumps;
2905:   PetscBool      isSeqSBAIJ;

2908:  #if defined(PETSC_USE_COMPLEX)
2909:   if (A->hermitian && !A->symmetric) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Hermitian CHOLESKY Factor is not supported");
2910:  #endif
2911:   MatCreate(PetscObjectComm((PetscObject)A),&B);
2912:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
2913:   PetscStrallocpy("mumps",&((PetscObject)B)->type_name);
2914:   MatSetUp(B);

2916:   PetscNewLog(B,&mumps);
2917:   PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);
2918:   if (isSeqSBAIJ) {
2919:     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
2920:   } else {
2921:     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
2922:   }

2924:   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2925:   B->ops->view                   = MatView_MUMPS;
2926:   B->ops->getinfo                = MatGetInfo_MUMPS;

2928:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);
2929:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);
2930:   PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);
2931:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);
2932:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);
2933:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);
2934:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);
2935:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);
2936:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);
2937:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);
2938:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);
2939:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverse_C",MatMumpsGetInverse_MUMPS);
2940:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverseTranspose_C",MatMumpsGetInverseTranspose_MUMPS);

2942:   B->factortype = MAT_FACTOR_CHOLESKY;
2943: #if defined(PETSC_USE_COMPLEX)
2944:   mumps->sym = 2;
2945: #else
2946:   if (A->spd_set && A->spd) mumps->sym = 1;
2947:   else                      mumps->sym = 2;
2948: #endif

2950:   /* set solvertype */
2951:   PetscFree(B->solvertype);
2952:   PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);

2954:   B->ops->destroy = MatDestroy_MUMPS;
2955:   B->data         = (void*)mumps;

2957:   PetscInitializeMUMPS(A,mumps);

2959:   *F = B;
2960:   return(0);
2961: }

2963: static PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
2964: {
2965:   Mat            B;
2967:   Mat_MUMPS      *mumps;
2968:   PetscBool      isSeqBAIJ;

2971:   /* Create the factorization matrix */
2972:   PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);
2973:   MatCreate(PetscObjectComm((PetscObject)A),&B);
2974:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
2975:   PetscStrallocpy("mumps",&((PetscObject)B)->type_name);
2976:   MatSetUp(B);

2978:   PetscNewLog(B,&mumps);
2979:   if (ftype == MAT_FACTOR_LU) {
2980:     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
2981:     B->factortype            = MAT_FACTOR_LU;
2982:     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
2983:     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
2984:     mumps->sym = 0;
2985:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n");

2987:   B->ops->view        = MatView_MUMPS;
2988:   B->ops->getinfo     = MatGetInfo_MUMPS;

2990:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);
2991:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);
2992:   PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);
2993:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);
2994:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);
2995:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);
2996:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);
2997:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);
2998:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);
2999:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);
3000:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);
3001:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverse_C",MatMumpsGetInverse_MUMPS);
3002:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverseTranspose_C",MatMumpsGetInverseTranspose_MUMPS);

3004:   /* set solvertype */
3005:   PetscFree(B->solvertype);
3006:   PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);

3008:   B->ops->destroy = MatDestroy_MUMPS;
3009:   B->data         = (void*)mumps;

3011:   PetscInitializeMUMPS(A,mumps);

3013:   *F = B;
3014:   return(0);
3015: }

3017: /* MatGetFactor for Seq and MPI SELL matrices */
3018: static PetscErrorCode MatGetFactor_sell_mumps(Mat A,MatFactorType ftype,Mat *F)
3019: {
3020:   Mat            B;
3022:   Mat_MUMPS      *mumps;
3023:   PetscBool      isSeqSELL;

3026:   /* Create the factorization matrix */
3027:   PetscObjectTypeCompare((PetscObject)A,MATSEQSELL,&isSeqSELL);
3028:   MatCreate(PetscObjectComm((PetscObject)A),&B);
3029:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
3030:   PetscStrallocpy("mumps",&((PetscObject)B)->type_name);
3031:   MatSetUp(B);

3033:   PetscNewLog(B,&mumps);

3035:   B->ops->view        = MatView_MUMPS;
3036:   B->ops->getinfo     = MatGetInfo_MUMPS;

3038:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);
3039:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);
3040:   PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);
3041:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);
3042:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);
3043:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);
3044:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);
3045:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);
3046:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);
3047:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);
3048:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);

3050:   if (ftype == MAT_FACTOR_LU) {
3051:     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
3052:     B->factortype            = MAT_FACTOR_LU;
3053:     if (isSeqSELL) mumps->ConvertToTriples = MatConvertToTriples_seqsell_seqaij;
3054:     else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"To be implemented");
3055:     mumps->sym = 0;
3056:   } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"To be implemented");

3058:   /* set solvertype */
3059:   PetscFree(B->solvertype);
3060:   PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);

3062:   B->ops->destroy = MatDestroy_MUMPS;
3063:   B->data         = (void*)mumps;

3065:   PetscInitializeMUMPS(A,mumps);

3067:   *F = B;
3068:   return(0);
3069: }

3071: PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_MUMPS(void)
3072: {

3076:   MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);
3077:   MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);
3078:   MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);
3079:   MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);
3080:   MatSolverTypeRegister(MATSOLVERMUMPS,MATMPISBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);
3081:   MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);
3082:   MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);
3083:   MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);
3084:   MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);
3085:   MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);
3086:   MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQSELL,MAT_FACTOR_LU,MatGetFactor_sell_mumps);
3087:   return(0);
3088: }