Actual source code: mumps.c

petsc-master 2017-05-24
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 <petscblaslapack.h>

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

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

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

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

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

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

100:   PetscErrorCode (*ConvertToTriples)(Mat, int, MatReuse, int*, int**, int**, PetscScalar**);
101: } Mat_MUMPS;

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

105: static PetscErrorCode MatMumpsResetSchur_Private(Mat_MUMPS* mumps)
106: {

110:   PetscFree2(mumps->id.listvar_schur,mumps->id.schur);
111:   PetscFree(mumps->id.redrhs);
112:   PetscFree(mumps->schur_sol);
113:   PetscFree(mumps->schur_pivots);
114:   PetscFree(mumps->schur_work);
115:   mumps->id.size_schur = 0;
116:   mumps->id.ICNTL(19) = 0;
117:   return(0);
118: }

120: static PetscErrorCode MatMumpsFactorSchur_Private(Mat_MUMPS* mumps)
121: {
122:   PetscBLASInt   B_N,B_ierr,B_slda;

126:   if (mumps->schur_factored) {
127:     return(0);
128:   }
129:   PetscBLASIntCast(mumps->id.size_schur,&B_N);
130:   PetscBLASIntCast(mumps->id.schur_lld,&B_slda);
131:   if (!mumps->sym) { /* MUMPS always return a full Schur matrix */
132:     if (!mumps->schur_pivots) {
133:       PetscMalloc1(B_N,&mumps->schur_pivots);
134:     }
135:     PetscFPTrapPush(PETSC_FP_TRAP_OFF);
136:     PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,(PetscScalar*)mumps->id.schur,&B_slda,mumps->schur_pivots,&B_ierr));
137:     PetscFPTrapPop();
138:     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
139:   } else { /* either full or lower-triangular (not packed) */
140:     char ord[2];
141:     if (mumps->id.ICNTL(19) == 2 || mumps->id.ICNTL(19) == 3) { /* lower triangular stored by columns or full matrix */
142:       sprintf(ord,"L");
143:     } else { /* ICNTL(19) == 1 lower triangular stored by rows */
144:       sprintf(ord,"U");
145:     }
146:     if (mumps->schur_sym == 2) {
147:       if (!mumps->schur_pivots) {
148:         PetscScalar  lwork;

150:         PetscMalloc1(B_N,&mumps->schur_pivots);
151:         mumps->schur_B_lwork=-1;
152:         PetscFPTrapPush(PETSC_FP_TRAP_OFF);
153:         PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_(ord,&B_N,(PetscScalar*)mumps->id.schur,&B_slda,mumps->schur_pivots,&lwork,&mumps->schur_B_lwork,&B_ierr));
154:         PetscFPTrapPop();
155:         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYTRF Lapack routine %d",(int)B_ierr);
156:         PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&mumps->schur_B_lwork);
157:         PetscMalloc1(mumps->schur_B_lwork,&mumps->schur_work);
158:       }
159:       PetscFPTrapPush(PETSC_FP_TRAP_OFF);
160:       PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_(ord,&B_N,(PetscScalar*)mumps->id.schur,&B_slda,mumps->schur_pivots,mumps->schur_work,&mumps->schur_B_lwork,&B_ierr));
161:       PetscFPTrapPop();
162:       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
163:     } else {
164:       PetscFPTrapPush(PETSC_FP_TRAP_OFF);
165:       PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_(ord,&B_N,(PetscScalar*)mumps->id.schur,&B_slda,&B_ierr));
166:       PetscFPTrapPop();
167:       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
168:     }
169:   }
170:   mumps->schur_factored = PETSC_TRUE;
171:   return(0);
172: }

174: static PetscErrorCode MatMumpsInvertSchur_Private(Mat_MUMPS* mumps)
175: {
176:   PetscBLASInt   B_N,B_ierr,B_slda;

180:   MatMumpsFactorSchur_Private(mumps);
181:   PetscBLASIntCast(mumps->id.size_schur,&B_N);
182:   PetscBLASIntCast(mumps->id.schur_lld,&B_slda);
183:   if (!mumps->sym) { /* MUMPS always return a full Schur matrix */
184:     if (!mumps->schur_work) {
185:       PetscScalar lwork;

187:       mumps->schur_B_lwork = -1;
188:       PetscFPTrapPush(PETSC_FP_TRAP_OFF);
189:       PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,(PetscScalar*)mumps->id.schur,&B_N,mumps->schur_pivots,&lwork,&mumps->schur_B_lwork,&B_ierr));
190:       PetscFPTrapPop();
191:       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
192:       PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&mumps->schur_B_lwork);
193:       PetscMalloc1(mumps->schur_B_lwork,&mumps->schur_work);
194:     }
195:     PetscFPTrapPush(PETSC_FP_TRAP_OFF);
196:     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,(PetscScalar*)mumps->id.schur,&B_N,mumps->schur_pivots,mumps->schur_work,&mumps->schur_B_lwork,&B_ierr));
197:     PetscFPTrapPop();
198:     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
199:   } else { /* either full or lower-triangular (not packed) */
200:     char ord[2];
201:     if (mumps->id.ICNTL(19) == 2 || mumps->id.ICNTL(19) == 3) { /* lower triangular stored by columns or full matrix */
202:       sprintf(ord,"L");
203:     } else { /* ICNTL(19) == 1 lower triangular stored by rows */
204:       sprintf(ord,"U");
205:     }
206:     if (mumps->schur_sym == 2) {
207:       PetscFPTrapPush(PETSC_FP_TRAP_OFF);
208:       PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_(ord,&B_N,(PetscScalar*)mumps->id.schur,&B_N,mumps->schur_pivots,mumps->schur_work,&B_ierr));
209:       PetscFPTrapPop();
210:       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
211:     } else {
212:       PetscFPTrapPush(PETSC_FP_TRAP_OFF);
213:       PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_(ord,&B_N,(PetscScalar*)mumps->id.schur,&B_N,&B_ierr));
214:       PetscFPTrapPop();
215:       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
216:     }
217:   }
218:   mumps->schur_inverted = PETSC_TRUE;
219:   return(0);
220: }

222: static PetscErrorCode MatMumpsSolveSchur_Private(Mat_MUMPS* mumps, PetscBool sol_in_redrhs)
223: {
224:   PetscBLASInt   B_N,B_Nrhs,B_ierr,B_slda,B_rlda;
225:   PetscScalar    one=1.,zero=0.;

229:   MatMumpsFactorSchur_Private(mumps);
230:   PetscBLASIntCast(mumps->id.size_schur,&B_N);
231:   PetscBLASIntCast(mumps->id.schur_lld,&B_slda);
232:   PetscBLASIntCast(mumps->id.nrhs,&B_Nrhs);
233:   PetscBLASIntCast(mumps->id.lredrhs,&B_rlda);
234:   if (mumps->schur_inverted) {
235:     PetscInt sizesol = B_Nrhs*B_N;
236:     if (!mumps->schur_sol || sizesol > mumps->schur_sizesol) {
237:       PetscFree(mumps->schur_sol);
238:       PetscMalloc1(sizesol,&mumps->schur_sol);
239:       mumps->schur_sizesol = sizesol;
240:     }
241:     if (!mumps->sym) {
242:       char type[2];
243:       if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
244:         if (!mumps->id.ICNTL(9)) { /* transpose solve */
245:           sprintf(type,"N");
246:         } else {
247:           sprintf(type,"T");
248:         }
249:       } else { /* stored by columns */
250:         if (!mumps->id.ICNTL(9)) { /* transpose solve */
251:           sprintf(type,"T");
252:         } else {
253:           sprintf(type,"N");
254:         }
255:       }
256:       PetscStackCallBLAS("BLASgemm",BLASgemm_(type,"N",&B_N,&B_Nrhs,&B_N,&one,(PetscScalar*)mumps->id.schur,&B_slda,(PetscScalar*)mumps->id.redrhs,&B_rlda,&zero,mumps->schur_sol,&B_rlda));
257:     } else {
258:       char ord[2];
259:       if (mumps->id.ICNTL(19) == 2 || mumps->id.ICNTL(19) == 3) { /* lower triangular stored by columns or full matrix */
260:         sprintf(ord,"L");
261:       } else { /* ICNTL(19) == 1 lower triangular stored by rows */
262:         sprintf(ord,"U");
263:       }
264:       PetscStackCallBLAS("BLASsymm",BLASsymm_("L",ord,&B_N,&B_Nrhs,&one,(PetscScalar*)mumps->id.schur,&B_slda,(PetscScalar*)mumps->id.redrhs,&B_rlda,&zero,mumps->schur_sol,&B_rlda));
265:     }
266:     if (sol_in_redrhs) {
267:       PetscMemcpy(mumps->id.redrhs,mumps->schur_sol,sizesol*sizeof(PetscScalar));
268:     }
269:   } else { /* Schur complement has not yet been inverted */
270:     MumpsScalar *orhs=NULL;

272:     if (!sol_in_redrhs) {
273:       PetscInt sizesol = B_Nrhs*B_N;
274:       if (!mumps->schur_sol || sizesol > mumps->schur_sizesol) {
275:         PetscFree(mumps->schur_sol);
276:         PetscMalloc1(sizesol,&mumps->schur_sol);
277:         mumps->schur_sizesol = sizesol;
278:       }
279:       orhs = mumps->id.redrhs;
280:       PetscMemcpy(mumps->schur_sol,mumps->id.redrhs,sizesol*sizeof(PetscScalar));
281:       mumps->id.redrhs = (MumpsScalar*)mumps->schur_sol;
282:     }
283:     if (!mumps->sym) { /* MUMPS always return a full Schur matrix */
284:       char type[2];
285:       if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
286:         if (!mumps->id.ICNTL(9)) { /* transpose solve */
287:           sprintf(type,"N");
288:         } else {
289:           sprintf(type,"T");
290:         }
291:       } else { /* stored by columns */
292:         if (!mumps->id.ICNTL(9)) { /* transpose solve */
293:           sprintf(type,"T");
294:         } else {
295:           sprintf(type,"N");
296:         }
297:       }
298:       PetscFPTrapPush(PETSC_FP_TRAP_OFF);
299:       PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_(type,&B_N,&B_Nrhs,(PetscScalar*)mumps->id.schur,&B_slda,mumps->schur_pivots,(PetscScalar*)mumps->id.redrhs,&B_rlda,&B_ierr));
300:       PetscFPTrapPop();
301:       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRS Lapack routine %d",(int)B_ierr);
302:     } else { /* either full or lower-triangular (not packed) */
303:       char ord[2];
304:       if (mumps->id.ICNTL(19) == 2 || mumps->id.ICNTL(19) == 3) { /* lower triangular stored by columns or full matrix */
305:         sprintf(ord,"L");
306:       } else { /* ICNTL(19) == 1 lower triangular stored by rows */
307:         sprintf(ord,"U");
308:       }
309:       if (mumps->schur_sym == 2) {
310:         PetscFPTrapPush(PETSC_FP_TRAP_OFF);
311:         PetscStackCallBLAS("LAPACKsytrs",LAPACKsytrs_(ord,&B_N,&B_Nrhs,(PetscScalar*)mumps->id.schur,&B_slda,mumps->schur_pivots,(PetscScalar*)mumps->id.redrhs,&B_rlda,&B_ierr));
312:         PetscFPTrapPop();
313:         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRS Lapack routine %d",(int)B_ierr);
314:       } else {
315:         PetscFPTrapPush(PETSC_FP_TRAP_OFF);
316:         PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_(ord,&B_N,&B_Nrhs,(PetscScalar*)mumps->id.schur,&B_slda,(PetscScalar*)mumps->id.redrhs,&B_rlda,&B_ierr));
317:         PetscFPTrapPop();
318:         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRS Lapack routine %d",(int)B_ierr);
319:       }
320:     }
321:     if (!sol_in_redrhs) {
322:       mumps->id.redrhs = orhs;
323:     }
324:   }
325:   return(0);
326: }

328: static PetscErrorCode MatMumpsHandleSchur_Private(Mat_MUMPS* mumps, PetscBool expansion)
329: {

333:   if (!mumps->id.ICNTL(19)) { /* do nothing when Schur complement has not been computed */
334:     return(0);
335:   }
336:   if (!expansion) { /* prepare for the condensation step */
337:     PetscInt sizeredrhs = mumps->id.nrhs*mumps->id.size_schur;
338:     /* allocate MUMPS internal array to store reduced right-hand sides */
339:     if (!mumps->id.redrhs || sizeredrhs > mumps->sizeredrhs) {
340:       PetscFree(mumps->id.redrhs);
341:       mumps->id.lredrhs = mumps->id.size_schur;
342:       PetscMalloc1(mumps->id.nrhs*mumps->id.lredrhs,&mumps->id.redrhs);
343:       mumps->sizeredrhs = mumps->id.nrhs*mumps->id.lredrhs;
344:     }
345:     mumps->id.ICNTL(26) = 1; /* condensation phase */
346:   } else { /* prepare for the expansion step */
347:     /* solve Schur complement (this has to be done by the MUMPS user, so basically us) */
348:     MatMumpsSolveSchur_Private(mumps,PETSC_TRUE);
349:     mumps->id.ICNTL(26) = 2; /* expansion phase */
350:     PetscMUMPS_c(&mumps->id);
351:     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));
352:     /* restore defaults */
353:     mumps->id.ICNTL(26) = -1;
354:     /* free MUMPS internal array for redrhs if we have solved for multiple rhs in order to save memory space */
355:     if (mumps->id.nrhs > 1) {
356:       PetscFree(mumps->id.redrhs);
357:       mumps->id.lredrhs = 0;
358:       mumps->sizeredrhs = 0;
359:     }
360:   }
361:   return(0);
362: }

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

367:   input:
368:     A       - matrix in aij,baij or sbaij (bs=1) format
369:     shift   - 0: C style output triple; 1: Fortran style output triple.
370:     reuse   - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple
371:               MAT_REUSE_MATRIX:   only the values in v array are updated
372:   output:
373:     nnz     - dim of r, c, and v (number of local nonzero entries of A)
374:     r, c, v - row and col index, matrix values (matrix triples)

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

380:  */

382: PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
383: {
384:   const PetscInt *ai,*aj,*ajj,M=A->rmap->n;
385:   PetscInt       nz,rnz,i,j;
387:   PetscInt       *row,*col;
388:   Mat_SeqAIJ     *aa=(Mat_SeqAIJ*)A->data;

391:   *v=aa->a;
392:   if (reuse == MAT_INITIAL_MATRIX) {
393:     nz   = aa->nz;
394:     ai   = aa->i;
395:     aj   = aa->j;
396:     *nnz = nz;
397:     PetscMalloc1(2*nz, &row);
398:     col  = row + nz;

400:     nz = 0;
401:     for (i=0; i<M; i++) {
402:       rnz = ai[i+1] - ai[i];
403:       ajj = aj + ai[i];
404:       for (j=0; j<rnz; j++) {
405:         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
406:       }
407:     }
408:     *r = row; *c = col;
409:   }
410:   return(0);
411: }

413: PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
414: {
415:   Mat_SeqBAIJ    *aa=(Mat_SeqBAIJ*)A->data;
416:   const PetscInt *ai,*aj,*ajj,bs2 = aa->bs2;
417:   PetscInt       bs,M,nz,idx=0,rnz,i,j,k,m;
419:   PetscInt       *row,*col;

422:   MatGetBlockSize(A,&bs);
423:   M = A->rmap->N/bs;
424:   *v = aa->a;
425:   if (reuse == MAT_INITIAL_MATRIX) {
426:     ai   = aa->i; aj = aa->j;
427:     nz   = bs2*aa->nz;
428:     *nnz = nz;
429:     PetscMalloc1(2*nz, &row);
430:     col  = row + nz;

432:     for (i=0; i<M; i++) {
433:       ajj = aj + ai[i];
434:       rnz = ai[i+1] - ai[i];
435:       for (k=0; k<rnz; k++) {
436:         for (j=0; j<bs; j++) {
437:           for (m=0; m<bs; m++) {
438:             row[idx]   = i*bs + m + shift;
439:             col[idx++] = bs*(ajj[k]) + j + shift;
440:           }
441:         }
442:       }
443:     }
444:     *r = row; *c = col;
445:   }
446:   return(0);
447: }

449: PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
450: {
451:   const PetscInt *ai, *aj,*ajj,M=A->rmap->n;
452:   PetscInt       nz,rnz,i,j;
454:   PetscInt       *row,*col;
455:   Mat_SeqSBAIJ   *aa=(Mat_SeqSBAIJ*)A->data;

458:   *v = aa->a;
459:   if (reuse == MAT_INITIAL_MATRIX) {
460:     nz   = aa->nz;
461:     ai   = aa->i;
462:     aj   = aa->j;
463:     *v   = aa->a;
464:     *nnz = nz;
465:     PetscMalloc1(2*nz, &row);
466:     col  = row + nz;

468:     nz = 0;
469:     for (i=0; i<M; i++) {
470:       rnz = ai[i+1] - ai[i];
471:       ajj = aj + ai[i];
472:       for (j=0; j<rnz; j++) {
473:         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
474:       }
475:     }
476:     *r = row; *c = col;
477:   }
478:   return(0);
479: }

481: PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
482: {
483:   const PetscInt    *ai,*aj,*ajj,*adiag,M=A->rmap->n;
484:   PetscInt          nz,rnz,i,j;
485:   const PetscScalar *av,*v1;
486:   PetscScalar       *val;
487:   PetscErrorCode    ierr;
488:   PetscInt          *row,*col;
489:   Mat_SeqAIJ        *aa=(Mat_SeqAIJ*)A->data;
490:   PetscBool         missing;

493:   ai   =aa->i; aj=aa->j;av=aa->a;
494:   adiag=aa->diag;
495:   MatMissingDiagonal_SeqAIJ(A,&missing,&i);
496:   if (reuse == MAT_INITIAL_MATRIX) {
497:     /* count nz in the uppper triangular part of A */
498:     nz = 0;
499:     if (missing) {
500:       for (i=0; i<M; i++) {
501:         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
502:           for (j=ai[i];j<ai[i+1];j++) {
503:             if (aj[j] < i) continue;
504:             nz++;
505:           }
506:         } else {
507:           nz += ai[i+1] - adiag[i];
508:         }
509:       }
510:     } else {
511:       for (i=0; i<M; i++) nz += ai[i+1] - adiag[i];
512:     }
513:     *nnz = nz;

515:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
516:     col  = row + nz;
517:     val  = (PetscScalar*)(col + nz);

519:     nz = 0;
520:     if (missing) {
521:       for (i=0; i<M; i++) {
522:         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
523:           for (j=ai[i];j<ai[i+1];j++) {
524:             if (aj[j] < i) continue;
525:             row[nz] = i+shift;
526:             col[nz] = aj[j]+shift;
527:             val[nz] = av[j];
528:             nz++;
529:           }
530:         } else {
531:           rnz = ai[i+1] - adiag[i];
532:           ajj = aj + adiag[i];
533:           v1  = av + adiag[i];
534:           for (j=0; j<rnz; j++) {
535:             row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
536:           }
537:         }
538:       }
539:     } else {
540:       for (i=0; i<M; i++) {
541:         rnz = ai[i+1] - adiag[i];
542:         ajj = aj + adiag[i];
543:         v1  = av + adiag[i];
544:         for (j=0; j<rnz; j++) {
545:           row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
546:         }
547:       }
548:     }
549:     *r = row; *c = col; *v = val;
550:   } else {
551:     nz = 0; val = *v;
552:     if (missing) {
553:       for (i=0; i <M; i++) {
554:         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
555:           for (j=ai[i];j<ai[i+1];j++) {
556:             if (aj[j] < i) continue;
557:             val[nz++] = av[j];
558:           }
559:         } else {
560:           rnz = ai[i+1] - adiag[i];
561:           v1  = av + adiag[i];
562:           for (j=0; j<rnz; j++) {
563:             val[nz++] = v1[j];
564:           }
565:         }
566:       }
567:     } else {
568:       for (i=0; i <M; i++) {
569:         rnz = ai[i+1] - adiag[i];
570:         v1  = av + adiag[i];
571:         for (j=0; j<rnz; j++) {
572:           val[nz++] = v1[j];
573:         }
574:       }
575:     }
576:   }
577:   return(0);
578: }

580: PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
581: {
582:   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
583:   PetscErrorCode    ierr;
584:   PetscInt          rstart,nz,i,j,jj,irow,countA,countB;
585:   PetscInt          *row,*col;
586:   const PetscScalar *av, *bv,*v1,*v2;
587:   PetscScalar       *val;
588:   Mat_MPISBAIJ      *mat = (Mat_MPISBAIJ*)A->data;
589:   Mat_SeqSBAIJ      *aa  = (Mat_SeqSBAIJ*)(mat->A)->data;
590:   Mat_SeqBAIJ       *bb  = (Mat_SeqBAIJ*)(mat->B)->data;

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

596:   garray = mat->garray;

598:   if (reuse == MAT_INITIAL_MATRIX) {
599:     nz   = aa->nz + bb->nz;
600:     *nnz = nz;
601:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
602:     col  = row + nz;
603:     val  = (PetscScalar*)(col + nz);

605:     *r = row; *c = col; *v = val;
606:   } else {
607:     row = *r; col = *c; val = *v;
608:   }

610:   jj = 0; irow = rstart;
611:   for (i=0; i<m; i++) {
612:     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
613:     countA = ai[i+1] - ai[i];
614:     countB = bi[i+1] - bi[i];
615:     bjj    = bj + bi[i];
616:     v1     = av + ai[i];
617:     v2     = bv + bi[i];

619:     /* A-part */
620:     for (j=0; j<countA; j++) {
621:       if (reuse == MAT_INITIAL_MATRIX) {
622:         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
623:       }
624:       val[jj++] = v1[j];
625:     }

627:     /* B-part */
628:     for (j=0; j < countB; j++) {
629:       if (reuse == MAT_INITIAL_MATRIX) {
630:         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
631:       }
632:       val[jj++] = v2[j];
633:     }
634:     irow++;
635:   }
636:   return(0);
637: }

639: PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
640: {
641:   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
642:   PetscErrorCode    ierr;
643:   PetscInt          rstart,nz,i,j,jj,irow,countA,countB;
644:   PetscInt          *row,*col;
645:   const PetscScalar *av, *bv,*v1,*v2;
646:   PetscScalar       *val;
647:   Mat_MPIAIJ        *mat = (Mat_MPIAIJ*)A->data;
648:   Mat_SeqAIJ        *aa  = (Mat_SeqAIJ*)(mat->A)->data;
649:   Mat_SeqAIJ        *bb  = (Mat_SeqAIJ*)(mat->B)->data;

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

655:   garray = mat->garray;

657:   if (reuse == MAT_INITIAL_MATRIX) {
658:     nz   = aa->nz + bb->nz;
659:     *nnz = nz;
660:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
661:     col  = row + nz;
662:     val  = (PetscScalar*)(col + nz);

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

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

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

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

698: PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
699: {
700:   Mat_MPIBAIJ       *mat    = (Mat_MPIBAIJ*)A->data;
701:   Mat_SeqBAIJ       *aa     = (Mat_SeqBAIJ*)(mat->A)->data;
702:   Mat_SeqBAIJ       *bb     = (Mat_SeqBAIJ*)(mat->B)->data;
703:   const PetscInt    *ai     = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj;
704:   const PetscInt    *garray = mat->garray,mbs=mat->mbs,rstart=A->rmap->rstart;
705:   const PetscInt    bs2=mat->bs2;
706:   PetscErrorCode    ierr;
707:   PetscInt          bs,nz,i,j,k,n,jj,irow,countA,countB,idx;
708:   PetscInt          *row,*col;
709:   const PetscScalar *av=aa->a, *bv=bb->a,*v1,*v2;
710:   PetscScalar       *val;

713:   MatGetBlockSize(A,&bs);
714:   if (reuse == MAT_INITIAL_MATRIX) {
715:     nz   = bs2*(aa->nz + bb->nz);
716:     *nnz = nz;
717:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
718:     col  = row + nz;
719:     val  = (PetscScalar*)(col + nz);

721:     *r = row; *c = col; *v = val;
722:   } else {
723:     row = *r; col = *c; val = *v;
724:   }

726:   jj = 0; irow = rstart;
727:   for (i=0; i<mbs; i++) {
728:     countA = ai[i+1] - ai[i];
729:     countB = bi[i+1] - bi[i];
730:     ajj    = aj + ai[i];
731:     bjj    = bj + bi[i];
732:     v1     = av + bs2*ai[i];
733:     v2     = bv + bs2*bi[i];

735:     idx = 0;
736:     /* A-part */
737:     for (k=0; k<countA; k++) {
738:       for (j=0; j<bs; j++) {
739:         for (n=0; n<bs; n++) {
740:           if (reuse == MAT_INITIAL_MATRIX) {
741:             row[jj] = irow + n + shift;
742:             col[jj] = rstart + bs*ajj[k] + j + shift;
743:           }
744:           val[jj++] = v1[idx++];
745:         }
746:       }
747:     }

749:     idx = 0;
750:     /* B-part */
751:     for (k=0; k<countB; k++) {
752:       for (j=0; j<bs; j++) {
753:         for (n=0; n<bs; n++) {
754:           if (reuse == MAT_INITIAL_MATRIX) {
755:             row[jj] = irow + n + shift;
756:             col[jj] = bs*garray[bjj[k]] + j + shift;
757:           }
758:           val[jj++] = v2[idx++];
759:         }
760:       }
761:     }
762:     irow += bs;
763:   }
764:   return(0);
765: }

767: PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
768: {
769:   const PetscInt    *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
770:   PetscErrorCode    ierr;
771:   PetscInt          rstart,nz,nza,nzb,i,j,jj,irow,countA,countB;
772:   PetscInt          *row,*col;
773:   const PetscScalar *av, *bv,*v1,*v2;
774:   PetscScalar       *val;
775:   Mat_MPIAIJ        *mat =  (Mat_MPIAIJ*)A->data;
776:   Mat_SeqAIJ        *aa  =(Mat_SeqAIJ*)(mat->A)->data;
777:   Mat_SeqAIJ        *bb  =(Mat_SeqAIJ*)(mat->B)->data;

780:   ai=aa->i; aj=aa->j; adiag=aa->diag;
781:   bi=bb->i; bj=bb->j; garray = mat->garray;
782:   av=aa->a; bv=bb->a;

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

786:   if (reuse == MAT_INITIAL_MATRIX) {
787:     nza = 0;    /* num of upper triangular entries in mat->A, including diagonals */
788:     nzb = 0;    /* num of upper triangular entries in mat->B */
789:     for (i=0; i<m; i++) {
790:       nza   += (ai[i+1] - adiag[i]);
791:       countB = bi[i+1] - bi[i];
792:       bjj    = bj + bi[i];
793:       for (j=0; j<countB; j++) {
794:         if (garray[bjj[j]] > rstart) nzb++;
795:       }
796:     }

798:     nz   = nza + nzb; /* total nz of upper triangular part of mat */
799:     *nnz = nz;
800:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
801:     col  = row + nz;
802:     val  = (PetscScalar*)(col + nz);

804:     *r = row; *c = col; *v = val;
805:   } else {
806:     row = *r; col = *c; val = *v;
807:   }

809:   jj = 0; irow = rstart;
810:   for (i=0; i<m; i++) {
811:     ajj    = aj + adiag[i];                 /* ptr to the beginning of the diagonal of this row */
812:     v1     = av + adiag[i];
813:     countA = ai[i+1] - adiag[i];
814:     countB = bi[i+1] - bi[i];
815:     bjj    = bj + bi[i];
816:     v2     = bv + bi[i];

818:     /* A-part */
819:     for (j=0; j<countA; j++) {
820:       if (reuse == MAT_INITIAL_MATRIX) {
821:         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
822:       }
823:       val[jj++] = v1[j];
824:     }

826:     /* B-part */
827:     for (j=0; j < countB; j++) {
828:       if (garray[bjj[j]] > rstart) {
829:         if (reuse == MAT_INITIAL_MATRIX) {
830:           row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
831:         }
832:         val[jj++] = v2[j];
833:       }
834:     }
835:     irow++;
836:   }
837:   return(0);
838: }

840: PetscErrorCode MatDestroy_MUMPS(Mat A)
841: {
842:   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->data;

846:   PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);
847:   VecScatterDestroy(&mumps->scat_rhs);
848:   VecScatterDestroy(&mumps->scat_sol);
849:   VecDestroy(&mumps->b_seq);
850:   VecDestroy(&mumps->x_seq);
851:   PetscFree(mumps->id.perm_in);
852:   PetscFree(mumps->irn);
853:   PetscFree(mumps->info);
854:   MatMumpsResetSchur_Private(mumps);
855:   mumps->id.job = JOB_END;
856:   PetscMUMPS_c(&mumps->id);
857:   MPI_Comm_free(&mumps->comm_mumps);
858:   PetscFree(A->data);

860:   /* clear composed functions */
861:   PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);
862:   PetscObjectComposeFunction((PetscObject)A,"MatFactorSetSchurIS_C",NULL);
863:   PetscObjectComposeFunction((PetscObject)A,"MatFactorInvertSchurComplement_C",NULL);
864:   PetscObjectComposeFunction((PetscObject)A,"MatFactorCreateSchurComplement_C",NULL);
865:   PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSchurComplement_C",NULL);
866:   PetscObjectComposeFunction((PetscObject)A,"MatFactorSolveSchurComplement_C",NULL);
867:   PetscObjectComposeFunction((PetscObject)A,"MatFactorSolveSchurComplementTranspose_C",NULL);
868:   PetscObjectComposeFunction((PetscObject)A,"MatFactorFactorizeSchurComplement_C",NULL);
869:   PetscObjectComposeFunction((PetscObject)A,"MatFactorSetSchurComplementSolverType_C",NULL);
870:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetIcntl_C",NULL);
871:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetIcntl_C",NULL);
872:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetCntl_C",NULL);
873:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetCntl_C",NULL);
874:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfo_C",NULL);
875:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfog_C",NULL);
876:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfo_C",NULL);
877:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfog_C",NULL);
878:   return(0);
879: }

881: PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x)
882: {
883:   Mat_MUMPS        *mumps=(Mat_MUMPS*)A->data;
884:   PetscScalar      *array;
885:   Vec              b_seq;
886:   IS               is_iden,is_petsc;
887:   PetscErrorCode   ierr;
888:   PetscInt         i;
889:   PetscBool        second_solve = PETSC_FALSE;
890:   static PetscBool cite1 = PETSC_FALSE,cite2 = PETSC_FALSE;

893:   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);
894:   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);

896:   if (A->factorerrortype) {
897:     PetscInfo2(A,"MatSolve is called with singular matrix factor, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
898:     VecSetInf(x);
899:     return(0);
900:   }

902:   mumps->id.nrhs = 1;
903:   b_seq          = mumps->b_seq;
904:   if (mumps->size > 1) {
905:     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
906:     VecScatterBegin(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);
907:     VecScatterEnd(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);
908:     if (!mumps->myid) {VecGetArray(b_seq,&array);}
909:   } else {  /* size == 1 */
910:     VecCopy(b,x);
911:     VecGetArray(x,&array);
912:   }
913:   if (!mumps->myid) { /* define rhs on the host */
914:     mumps->id.nrhs = 1;
915:     mumps->id.rhs = (MumpsScalar*)array;
916:   }

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

936:   /* handle expansion step of Schur complement (if any) */
937:   if (second_solve) {
938:     MatMumpsHandleSchur_Private(mumps,PETSC_TRUE);
939:   }

941:   if (mumps->size > 1) { /* convert mumps distributed solution to petsc mpi x */
942:     if (mumps->scat_sol && mumps->ICNTL9_pre != mumps->id.ICNTL(9)) {
943:       /* when id.ICNTL(9) changes, the contents of lsol_loc may change (not its size, lsol_loc), recreates scat_sol */
944:       VecScatterDestroy(&mumps->scat_sol);
945:     }
946:     if (!mumps->scat_sol) { /* create scatter scat_sol */
947:       ISCreateStride(PETSC_COMM_SELF,mumps->id.lsol_loc,0,1,&is_iden); /* from */
948:       for (i=0; i<mumps->id.lsol_loc; i++) {
949:         mumps->id.isol_loc[i] -= 1; /* change Fortran style to C style */
950:       }
951:       ISCreateGeneral(PETSC_COMM_SELF,mumps->id.lsol_loc,mumps->id.isol_loc,PETSC_COPY_VALUES,&is_petsc);  /* to */
952:       VecScatterCreate(mumps->x_seq,is_iden,x,is_petsc,&mumps->scat_sol);
953:       ISDestroy(&is_iden);
954:       ISDestroy(&is_petsc);

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

959:     VecScatterBegin(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);
960:     VecScatterEnd(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);
961:   }
962:   return(0);
963: }

965: PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x)
966: {
967:   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->data;

971:   mumps->id.ICNTL(9) = 0;
972:   MatSolve_MUMPS(A,b,x);
973:   mumps->id.ICNTL(9) = 1;
974:   return(0);
975: }

977: PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X)
978: {
980:   Mat            Bt = NULL;
981:   PetscBool      flg, flgT;
982:   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->data;
983:   PetscInt       i,nrhs,M;
984:   PetscScalar    *array,*bray;

987:   PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
988:   PetscObjectTypeCompare((PetscObject)B,MATTRANSPOSEMAT,&flgT);
989:   if (flgT) {
990:     if (mumps->size > 1) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
991:     MatTransposeGetMat(B,&Bt);
992:   } else {
993:     if (!flg) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
994:   }
995:   PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
996:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
997:   if (B->rmap->n != X->rmap->n) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B and X must have same row distribution");

999:   MatGetSize(B,&M,&nrhs);
1000:   mumps->id.nrhs = nrhs;
1001:   mumps->id.lrhs = M;

1003:   if (mumps->size == 1) {
1004:     PetscScalar *aa;
1005:     PetscInt    spnr,*ia,*ja;
1006:     PetscBool   second_solve = PETSC_FALSE;

1008:     /* copy B to X */
1009:     MatDenseGetArray(X,&array);
1010:     mumps->id.rhs = (MumpsScalar*)array;
1011:     if (!Bt) {
1012:       MatDenseGetArray(B,&bray);
1013:       PetscMemcpy(array,bray,M*nrhs*sizeof(PetscScalar));
1014:       MatDenseRestoreArray(B,&bray);
1015:     } else {
1016:       PetscBool done;

1018:       MatSeqAIJGetArray(Bt,&aa);
1019:       MatGetRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&done);
1020:       if (!done) SETERRQ(PetscObjectComm((PetscObject)Bt),PETSC_ERR_ARG_WRONG,"Cannot get IJ structure");
1021:       mumps->id.irhs_ptr    = ia;
1022:       mumps->id.irhs_sparse = ja;
1023:       mumps->id.nz_rhs      = ia[spnr] - 1;
1024:       mumps->id.rhs_sparse  = (MumpsScalar*)aa;
1025:       mumps->id.ICNTL(20)   = 1;
1026:     }
1027:     /* handle condensation step of Schur complement (if any) */
1028:     if (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2) {
1029:       second_solve = PETSC_TRUE;
1030:       MatMumpsHandleSchur_Private(mumps,PETSC_FALSE);
1031:     }
1032:     /* solve phase */
1033:     /*-------------*/
1034:     mumps->id.job = JOB_SOLVE;
1035:     PetscMUMPS_c(&mumps->id);
1036:     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));

1038:     /* handle expansion step of Schur complement (if any) */
1039:     if (second_solve) {
1040:       MatMumpsHandleSchur_Private(mumps,PETSC_TRUE);
1041:     }
1042:     if (Bt) {
1043:       PetscBool done;

1045:       MatSeqAIJRestoreArray(Bt,&aa);
1046:       MatRestoreRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&done);
1047:       if (!done) SETERRQ(PetscObjectComm((PetscObject)Bt),PETSC_ERR_ARG_WRONG,"Cannot restore IJ structure");
1048:       mumps->id.ICNTL(20) = 0;
1049:     }
1050:     MatDenseRestoreArray(X,&array);
1051:   } else {  /*--------- parallel case --------*/
1052:     PetscInt       lsol_loc,nlsol_loc,*isol_loc,*idx,*iidx,*idxx,*isol_loc_save;
1053:     MumpsScalar    *sol_loc,*sol_loc_save;
1054:     IS             is_to,is_from;
1055:     PetscInt       k,proc,j,m;
1056:     const PetscInt *rstart;
1057:     Vec            v_mpi,b_seq,x_seq;
1058:     VecScatter     scat_rhs,scat_sol;

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

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

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

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

1074:     /* copy rhs matrix B into vector v_mpi */
1075:     MatGetLocalSize(B,&m,NULL);
1076:     MatDenseGetArray(B,&bray);
1077:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)B),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi);
1078:     MatDenseRestoreArray(B,&bray);

1080:     /* scatter v_mpi to b_seq because MUMPS only supports centralized rhs */
1081:     /* idx: maps from k-th index of v_mpi to (i,j)-th global entry of B;
1082:       iidx: inverse of idx, will be used by scattering xx_seq -> X       */
1083:     PetscMalloc2(nrhs*M,&idx,nrhs*M,&iidx);
1084:     MatGetOwnershipRanges(B,&rstart);
1085:     k = 0;
1086:     for (proc=0; proc<mumps->size; proc++){
1087:       for (j=0; j<nrhs; j++){
1088:         for (i=rstart[proc]; i<rstart[proc+1]; i++){
1089:           iidx[j*M + i] = k;
1090:           idx[k++]      = j*M + i;
1091:         }
1092:       }
1093:     }

1095:     if (!mumps->myid) {
1096:       VecCreateSeq(PETSC_COMM_SELF,nrhs*M,&b_seq);
1097:       ISCreateGeneral(PETSC_COMM_SELF,nrhs*M,idx,PETSC_COPY_VALUES,&is_to);
1098:       ISCreateStride(PETSC_COMM_SELF,nrhs*M,0,1,&is_from);
1099:     } else {
1100:       VecCreateSeq(PETSC_COMM_SELF,0,&b_seq);
1101:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_to);
1102:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_from);
1103:     }
1104:     VecScatterCreate(v_mpi,is_from,b_seq,is_to,&scat_rhs);
1105:     VecScatterBegin(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);
1106:     ISDestroy(&is_to);
1107:     ISDestroy(&is_from);
1108:     VecScatterEnd(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);

1110:     if (!mumps->myid) { /* define rhs on the host */
1111:       VecGetArray(b_seq,&bray);
1112:       mumps->id.rhs = (MumpsScalar*)bray;
1113:       VecRestoreArray(b_seq,&bray);
1114:     }

1116:     /* solve phase */
1117:     /*-------------*/
1118:     mumps->id.job = JOB_SOLVE;
1119:     PetscMUMPS_c(&mumps->id);
1120:     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));

1122:     /* scatter mumps distributed solution to petsc vector v_mpi, which shares local arrays with solution matrix X */
1123:     MatDenseGetArray(X,&array);
1124:     VecPlaceArray(v_mpi,array);
1125: 
1126:     /* create scatter scat_sol */
1127:     PetscMalloc1(nlsol_loc,&idxx);
1128:     ISCreateStride(PETSC_COMM_SELF,nlsol_loc,0,1,&is_from);
1129:     for (i=0; i<lsol_loc; i++) {
1130:       isol_loc[i] -= 1; /* change Fortran style to C style */
1131:       idxx[i] = iidx[isol_loc[i]];
1132:       for (j=1; j<nrhs; j++){
1133:         idxx[j*lsol_loc+i] = iidx[isol_loc[i]+j*M];
1134:       }
1135:     }
1136:     ISCreateGeneral(PETSC_COMM_SELF,nlsol_loc,idxx,PETSC_COPY_VALUES,&is_to);
1137:     VecScatterCreate(x_seq,is_from,v_mpi,is_to,&scat_sol);
1138:     VecScatterBegin(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);
1139:     ISDestroy(&is_from);
1140:     ISDestroy(&is_to);
1141:     VecScatterEnd(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);
1142:     MatDenseRestoreArray(X,&array);

1144:     /* free spaces */
1145:     mumps->id.sol_loc = sol_loc_save;
1146:     mumps->id.isol_loc = isol_loc_save;

1148:     PetscFree2(sol_loc,isol_loc);
1149:     PetscFree2(idx,iidx);
1150:     PetscFree(idxx);
1151:     VecDestroy(&x_seq);
1152:     VecDestroy(&v_mpi);
1153:     VecDestroy(&b_seq);
1154:     VecScatterDestroy(&scat_rhs);
1155:     VecScatterDestroy(&scat_sol);
1156:   }
1157:   return(0);
1158: }

1160: #if !defined(PETSC_USE_COMPLEX)
1161: /*
1162:   input:
1163:    F:        numeric factor
1164:   output:
1165:    nneg:     total number of negative pivots
1166:    nzero:    total number of zero pivots
1167:    npos:     (global dimension of F) - nneg - nzero
1168: */
1169: PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
1170: {
1171:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1173:   PetscMPIInt    size;

1176:   MPI_Comm_size(PetscObjectComm((PetscObject)F),&size);
1177:   /* 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 */
1178:   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));

1180:   if (nneg) *nneg = mumps->id.INFOG(12);
1181:   if (nzero || npos) {
1182:     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");
1183:     if (nzero) *nzero = mumps->id.INFOG(28);
1184:     if (npos) *npos   = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28));
1185:   }
1186:   return(0);
1187: }
1188: #endif

1190: PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
1191: {
1192:   Mat_MUMPS      *mumps =(Mat_MUMPS*)(F)->data;
1194:   PetscBool      isMPIAIJ;

1197:   if (mumps->id.INFOG(1) < 0) {
1198:     if (mumps->id.INFOG(1) == -6) {
1199:       PetscInfo2(A,"MatFactorNumeric is called with singular matrix structure, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1200:     }
1201:     PetscInfo2(A,"MatFactorNumeric is called after analysis phase fails, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1202:     return(0);
1203:   }

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

1207:   /* numerical factorization phase */
1208:   /*-------------------------------*/
1209:   mumps->id.job = JOB_FACTNUMERIC;
1210:   if (!mumps->id.ICNTL(18)) { /* A is centralized */
1211:     if (!mumps->myid) {
1212:       mumps->id.a = (MumpsScalar*)mumps->val;
1213:     }
1214:   } else {
1215:     mumps->id.a_loc = (MumpsScalar*)mumps->val;
1216:   }
1217:   PetscMUMPS_c(&mumps->id);
1218:   if (mumps->id.INFOG(1) < 0) {
1219:     if (A->erroriffailure) {
1220:       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));
1221:     } else {
1222:       if (mumps->id.INFOG(1) == -10) { /* numerically singular matrix */
1223:         PetscInfo2(F,"matrix is numerically singular, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1224:         F->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1225:       } else if (mumps->id.INFOG(1) == -13) {
1226:         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));
1227:         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1228:       } else if (mumps->id.INFOG(1) == -8 || mumps->id.INFOG(1) == -9 || (-16 < mumps->id.INFOG(1) && mumps->id.INFOG(1) < -10) ) {
1229:         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));
1230:         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1231:       } else {
1232:         PetscInfo2(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1233:         F->factorerrortype = MAT_FACTOR_OTHER;
1234:       }
1235:     }
1236:   }
1237:   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));

1239:   (F)->assembled        = PETSC_TRUE;
1240:   mumps->matstruc       = SAME_NONZERO_PATTERN;
1241:   mumps->schur_factored = PETSC_FALSE;
1242:   mumps->schur_inverted = PETSC_FALSE;

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

1247:   if (mumps->size > 1) {
1248:     PetscInt    lsol_loc;
1249:     PetscScalar *sol_loc;

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

1253:     /* distributed solution; Create x_seq=sol_loc for repeated use */
1254:     if (mumps->x_seq) {
1255:       VecScatterDestroy(&mumps->scat_sol);
1256:       PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);
1257:       VecDestroy(&mumps->x_seq);
1258:     }
1259:     lsol_loc = mumps->id.INFO(23); /* length of sol_loc */
1260:     PetscMalloc2(lsol_loc,&sol_loc,lsol_loc,&mumps->id.isol_loc);
1261:     mumps->id.lsol_loc = lsol_loc;
1262:     mumps->id.sol_loc = (MumpsScalar*)sol_loc;
1263:     VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq);
1264:   }
1265:   return(0);
1266: }

1268: /* Sets MUMPS options from the options database */
1269: PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A)
1270: {
1271:   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1273:   PetscInt       icntl,info[40],i,ninfo=40;
1274:   PetscBool      flg;

1277:   PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MUMPS Options","Mat");
1278:   PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);
1279:   if (flg) mumps->id.ICNTL(1) = icntl;
1280:   PetscOptionsInt("-mat_mumps_icntl_2","ICNTL(2): output stream for diagnostic printing, statistics, and warning","None",mumps->id.ICNTL(2),&icntl,&flg);
1281:   if (flg) mumps->id.ICNTL(2) = icntl;
1282:   PetscOptionsInt("-mat_mumps_icntl_3","ICNTL(3): output stream for global information, collected on the host","None",mumps->id.ICNTL(3),&icntl,&flg);
1283:   if (flg) mumps->id.ICNTL(3) = icntl;

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

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

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

1298:   PetscOptionsInt("-mat_mumps_icntl_8","ICNTL(8): scaling strategy (-2 to 8 or 77)","None",mumps->id.ICNTL(8),&mumps->id.ICNTL(8),NULL);
1299:   /* 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() */
1300:   PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",mumps->id.ICNTL(10),&mumps->id.ICNTL(10),NULL);
1301:   PetscOptionsInt("-mat_mumps_icntl_11","ICNTL(11): statistics related to an error analysis (via -ksp_view)","None",mumps->id.ICNTL(11),&mumps->id.ICNTL(11),NULL);
1302:   PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): an ordering strategy for symmetric matrices (0 to 3)","None",mumps->id.ICNTL(12),&mumps->id.ICNTL(12),NULL);
1303:   PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting","None",mumps->id.ICNTL(13),&mumps->id.ICNTL(13),NULL);
1304:   PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage increase in the estimated working space","None",mumps->id.ICNTL(14),&mumps->id.ICNTL(14),NULL);
1305:   PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): computes the Schur complement","None",mumps->id.ICNTL(19),&mumps->id.ICNTL(19),NULL);
1306:   if (mumps->id.ICNTL(19) <= 0 || mumps->id.ICNTL(19) > 3) { /* reset any schur data (if any) */
1307:     MatMumpsResetSchur_Private(mumps);
1308:   }
1309:   /* PetscOptionsInt("-mat_mumps_icntl_20","ICNTL(20): the format (dense or sparse) of the right-hand sides","None",mumps->id.ICNTL(20),&mumps->id.ICNTL(20),NULL); -- sparse rhs is not supported in PETSc API */
1310:   /* PetscOptionsInt("-mat_mumps_icntl_21","ICNTL(21): the distribution (centralized or distributed) of the solution vectors","None",mumps->id.ICNTL(21),&mumps->id.ICNTL(21),NULL); we only use distributed solution vector */

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

1319:   PetscOptionsInt("-mat_mumps_icntl_25","ICNTL(25): compute a solution of a deficient matrix and a null space basis","None",mumps->id.ICNTL(25),&mumps->id.ICNTL(25),NULL);
1320:   PetscOptionsInt("-mat_mumps_icntl_26","ICNTL(26): drives the solution phase if a Schur complement matrix","None",mumps->id.ICNTL(26),&mumps->id.ICNTL(26),NULL);
1321:   PetscOptionsInt("-mat_mumps_icntl_27","ICNTL(27): the blocking size for multiple right-hand sides","None",mumps->id.ICNTL(27),&mumps->id.ICNTL(27),NULL);
1322:   PetscOptionsInt("-mat_mumps_icntl_28","ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering","None",mumps->id.ICNTL(28),&mumps->id.ICNTL(28),NULL);
1323:   PetscOptionsInt("-mat_mumps_icntl_29","ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis","None",mumps->id.ICNTL(29),&mumps->id.ICNTL(29),NULL);
1324:   PetscOptionsInt("-mat_mumps_icntl_30","ICNTL(30): compute user-specified set of entries in inv(A)","None",mumps->id.ICNTL(30),&mumps->id.ICNTL(30),NULL);
1325:   PetscOptionsInt("-mat_mumps_icntl_31","ICNTL(31): indicates which factors may be discarded during factorization","None",mumps->id.ICNTL(31),&mumps->id.ICNTL(31),NULL);
1326:   /* PetscOptionsInt("-mat_mumps_icntl_32","ICNTL(32): performs the forward elemination of the right-hand sides during factorization","None",mumps->id.ICNTL(32),&mumps->id.ICNTL(32),NULL);  -- not supported by PETSc API */
1327:   PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),NULL);

1329:   PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),NULL);
1330:   PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",mumps->id.CNTL(2),&mumps->id.CNTL(2),NULL);
1331:   PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),NULL);
1332:   PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",mumps->id.CNTL(4),&mumps->id.CNTL(4),NULL);
1333:   PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",mumps->id.CNTL(5),&mumps->id.CNTL(5),NULL);

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

1337:   PetscOptionsIntArray("-mat_mumps_view_info","request INFO local to each processor","",info,&ninfo,NULL);
1338:   if (ninfo) {
1339:     if (ninfo > 40) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"number of INFO %d must <= 40\n",ninfo);
1340:     PetscMalloc1(ninfo,&mumps->info);
1341:     mumps->ninfo = ninfo;
1342:     for (i=0; i<ninfo; i++) {
1343:       if (info[i] < 0 || info[i]>40) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index of INFO %d must between 1 and 40\n",ninfo);
1344:       else  mumps->info[i] = info[i];
1345:     }
1346:   }

1348:   PetscOptionsEnd();
1349:   return(0);
1350: }

1352: PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps)
1353: {

1357:   MPI_Comm_rank(PetscObjectComm((PetscObject)A), &mumps->myid);
1358:   MPI_Comm_size(PetscObjectComm((PetscObject)A),&mumps->size);
1359:   MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(mumps->comm_mumps));

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

1363:   mumps->id.job = JOB_INIT;
1364:   mumps->id.par = 1;  /* host participates factorizaton and solve */
1365:   mumps->id.sym = mumps->sym;
1366:   PetscMUMPS_c(&mumps->id);

1368:   mumps->scat_rhs     = NULL;
1369:   mumps->scat_sol     = NULL;

1371:   /* set PETSc-MUMPS default options - override MUMPS default */
1372:   mumps->id.ICNTL(3) = 0;
1373:   mumps->id.ICNTL(4) = 0;
1374:   if (mumps->size == 1) {
1375:     mumps->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
1376:   } else {
1377:     mumps->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
1378:     mumps->id.ICNTL(20) = 0;   /* rhs is in dense format */
1379:     mumps->id.ICNTL(21) = 1;   /* distributed solution */
1380:   }

1382:   /* schur */
1383:   mumps->id.size_schur      = 0;
1384:   mumps->id.listvar_schur   = NULL;
1385:   mumps->id.schur           = NULL;
1386:   mumps->sizeredrhs         = 0;
1387:   mumps->schur_pivots       = NULL;
1388:   mumps->schur_work         = NULL;
1389:   mumps->schur_sol          = NULL;
1390:   mumps->schur_sizesol      = 0;
1391:   mumps->schur_factored     = PETSC_FALSE;
1392:   mumps->schur_inverted     = PETSC_FALSE;
1393:   mumps->schur_sym          = mumps->id.sym;
1394:   return(0);
1395: }

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

1402:   if (mumps->id.INFOG(1) < 0) {
1403:     if (A->erroriffailure) {
1404:       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
1405:     } else {
1406:       if (mumps->id.INFOG(1) == -6) {
1407:         PetscInfo2(F,"matrix is singular in structure, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1408:         F->factorerrortype = MAT_FACTOR_STRUCT_ZEROPIVOT;
1409:       } else if (mumps->id.INFOG(1) == -5 || mumps->id.INFOG(1) == -7) {
1410:         PetscInfo2(F,"problem of workspace, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1411:         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1412:       } else {
1413:         PetscInfo2(F,"Error reported by MUMPS in analysis phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1414:         F->factorerrortype = MAT_FACTOR_OTHER;
1415:       }
1416:     }
1417:   }
1418:   return(0);
1419: }

1421: /* Note Petsc r(=c) permutation is used when mumps->id.ICNTL(7)==1 with centralized assembled matrix input; otherwise r and c are ignored */
1422: PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1423: {
1424:   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1426:   Vec            b;
1427:   IS             is_iden;
1428:   const PetscInt M = A->rmap->N;

1431:   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;

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

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

1438:   /* analysis phase */
1439:   /*----------------*/
1440:   mumps->id.job = JOB_FACTSYMBOLIC;
1441:   mumps->id.n   = M;
1442:   switch (mumps->id.ICNTL(18)) {
1443:   case 0:  /* centralized assembled matrix input */
1444:     if (!mumps->myid) {
1445:       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1446:       if (mumps->id.ICNTL(6)>1) {
1447:         mumps->id.a = (MumpsScalar*)mumps->val;
1448:       }
1449:       if (mumps->id.ICNTL(7) == 1) { /* use user-provide matrix ordering - assuming r = c ordering */
1450:         /*
1451:         PetscBool      flag;
1452:         ISEqual(r,c,&flag);
1453:         if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"row_perm != col_perm");
1454:         ISView(r,PETSC_VIEWER_STDOUT_SELF);
1455:          */
1456:         if (!mumps->myid) {
1457:           const PetscInt *idx;
1458:           PetscInt       i,*perm_in;

1460:           PetscMalloc1(M,&perm_in);
1461:           ISGetIndices(r,&idx);

1463:           mumps->id.perm_in = perm_in;
1464:           for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */
1465:           ISRestoreIndices(r,&idx);
1466:         }
1467:       }
1468:     }
1469:     break;
1470:   case 3:  /* distributed assembled matrix input (size>1) */
1471:     mumps->id.nz_loc = mumps->nz;
1472:     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1473:     if (mumps->id.ICNTL(6)>1) {
1474:       mumps->id.a_loc = (MumpsScalar*)mumps->val;
1475:     }
1476:     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1477:     if (!mumps->myid) {
1478:       VecCreateSeq(PETSC_COMM_SELF,A->rmap->N,&mumps->b_seq);
1479:       ISCreateStride(PETSC_COMM_SELF,A->rmap->N,0,1,&is_iden);
1480:     } else {
1481:       VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);
1482:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);
1483:     }
1484:     MatCreateVecs(A,NULL,&b);
1485:     VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);
1486:     ISDestroy(&is_iden);
1487:     VecDestroy(&b);
1488:     break;
1489:   }
1490:   PetscMUMPS_c(&mumps->id);
1491:   MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);

1493:   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1494:   F->ops->solve           = MatSolve_MUMPS;
1495:   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1496:   F->ops->matsolve        = MatMatSolve_MUMPS;
1497:   return(0);
1498: }

1500: /* Note the Petsc r and c permutations are ignored */
1501: PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1502: {
1503:   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1505:   Vec            b;
1506:   IS             is_iden;
1507:   const PetscInt M = A->rmap->N;

1510:   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;

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

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

1517:   /* analysis phase */
1518:   /*----------------*/
1519:   mumps->id.job = JOB_FACTSYMBOLIC;
1520:   mumps->id.n   = M;
1521:   switch (mumps->id.ICNTL(18)) {
1522:   case 0:  /* centralized assembled matrix input */
1523:     if (!mumps->myid) {
1524:       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1525:       if (mumps->id.ICNTL(6)>1) {
1526:         mumps->id.a = (MumpsScalar*)mumps->val;
1527:       }
1528:     }
1529:     break;
1530:   case 3:  /* distributed assembled matrix input (size>1) */
1531:     mumps->id.nz_loc = mumps->nz;
1532:     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1533:     if (mumps->id.ICNTL(6)>1) {
1534:       mumps->id.a_loc = (MumpsScalar*)mumps->val;
1535:     }
1536:     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1537:     if (!mumps->myid) {
1538:       VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);
1539:       ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);
1540:     } else {
1541:       VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);
1542:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);
1543:     }
1544:     MatCreateVecs(A,NULL,&b);
1545:     VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);
1546:     ISDestroy(&is_iden);
1547:     VecDestroy(&b);
1548:     break;
1549:   }
1550:   PetscMUMPS_c(&mumps->id);
1551:   MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);

1553:   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1554:   F->ops->solve           = MatSolve_MUMPS;
1555:   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1556:   return(0);
1557: }

1559: /* Note the Petsc r permutation and factor info are ignored */
1560: PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
1561: {
1562:   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1564:   Vec            b;
1565:   IS             is_iden;
1566:   const PetscInt M = A->rmap->N;

1569:   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;

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

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

1576:   /* analysis phase */
1577:   /*----------------*/
1578:   mumps->id.job = JOB_FACTSYMBOLIC;
1579:   mumps->id.n   = M;
1580:   switch (mumps->id.ICNTL(18)) {
1581:   case 0:  /* centralized assembled matrix input */
1582:     if (!mumps->myid) {
1583:       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1584:       if (mumps->id.ICNTL(6)>1) {
1585:         mumps->id.a = (MumpsScalar*)mumps->val;
1586:       }
1587:     }
1588:     break;
1589:   case 3:  /* distributed assembled matrix input (size>1) */
1590:     mumps->id.nz_loc = mumps->nz;
1591:     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1592:     if (mumps->id.ICNTL(6)>1) {
1593:       mumps->id.a_loc = (MumpsScalar*)mumps->val;
1594:     }
1595:     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1596:     if (!mumps->myid) {
1597:       VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);
1598:       ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);
1599:     } else {
1600:       VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);
1601:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);
1602:     }
1603:     MatCreateVecs(A,NULL,&b);
1604:     VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);
1605:     ISDestroy(&is_iden);
1606:     VecDestroy(&b);
1607:     break;
1608:   }
1609:   PetscMUMPS_c(&mumps->id);
1610:   MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);

1612:   F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
1613:   F->ops->solve                 = MatSolve_MUMPS;
1614:   F->ops->solvetranspose        = MatSolve_MUMPS;
1615:   F->ops->matsolve              = MatMatSolve_MUMPS;
1616: #if defined(PETSC_USE_COMPLEX)
1617:   F->ops->getinertia = NULL;
1618: #else
1619:   F->ops->getinertia = MatGetInertia_SBAIJMUMPS;
1620: #endif
1621:   return(0);
1622: }

1624: PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
1625: {
1626:   PetscErrorCode    ierr;
1627:   PetscBool         iascii;
1628:   PetscViewerFormat format;
1629:   Mat_MUMPS         *mumps=(Mat_MUMPS*)A->data;

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

1635:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1636:   if (iascii) {
1637:     PetscViewerGetFormat(viewer,&format);
1638:     if (format == PETSC_VIEWER_ASCII_INFO) {
1639:       PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");
1640:       PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                   %d \n",mumps->id.sym);
1641:       PetscViewerASCIIPrintf(viewer,"  PAR (host participation):            %d \n",mumps->id.par);
1642:       PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):         %d \n",mumps->id.ICNTL(1));
1643:       PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg): %d \n",mumps->id.ICNTL(2));
1644:       PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):   %d \n",mumps->id.ICNTL(3));
1645:       PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):        %d \n",mumps->id.ICNTL(4));
1646:       PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):         %d \n",mumps->id.ICNTL(5));
1647:       PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):        %d \n",mumps->id.ICNTL(6));
1648:       PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (sequential matrix ordering):%d \n",mumps->id.ICNTL(7));
1649:       PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scaling strategy):        %d \n",mumps->id.ICNTL(8));
1650:       PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements):  %d \n",mumps->id.ICNTL(10));
1651:       PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):          %d \n",mumps->id.ICNTL(11));
1652:       if (mumps->id.ICNTL(11)>0) {
1653:         PetscViewerASCIIPrintf(viewer,"    RINFOG(4) (inf norm of input mat):        %g\n",mumps->id.RINFOG(4));
1654:         PetscViewerASCIIPrintf(viewer,"    RINFOG(5) (inf norm of solution):         %g\n",mumps->id.RINFOG(5));
1655:         PetscViewerASCIIPrintf(viewer,"    RINFOG(6) (inf norm of residual):         %g\n",mumps->id.RINFOG(6));
1656:         PetscViewerASCIIPrintf(viewer,"    RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8));
1657:         PetscViewerASCIIPrintf(viewer,"    RINFOG(9) (error estimate):               %g \n",mumps->id.RINFOG(9));
1658:         PetscViewerASCIIPrintf(viewer,"    RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11));
1659:       }
1660:       PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",mumps->id.ICNTL(12));
1661:       PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",mumps->id.ICNTL(13));
1662:       PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",mumps->id.ICNTL(14));
1663:       /* ICNTL(15-17) not used */
1664:       PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",mumps->id.ICNTL(18));
1665:       PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Schur complement info):                       %d \n",mumps->id.ICNTL(19));
1666:       PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",mumps->id.ICNTL(20));
1667:       PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d \n",mumps->id.ICNTL(21));
1668:       PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",mumps->id.ICNTL(22));
1669:       PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",mumps->id.ICNTL(23));

1671:       PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",mumps->id.ICNTL(24));
1672:       PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",mumps->id.ICNTL(25));
1673:       PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",mumps->id.ICNTL(26));
1674:       PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",mumps->id.ICNTL(27));
1675:       PetscViewerASCIIPrintf(viewer,"  ICNTL(28) (use parallel or sequential ordering):        %d \n",mumps->id.ICNTL(28));
1676:       PetscViewerASCIIPrintf(viewer,"  ICNTL(29) (parallel ordering):                          %d \n",mumps->id.ICNTL(29));

1678:       PetscViewerASCIIPrintf(viewer,"  ICNTL(30) (user-specified set of entries in inv(A)):    %d \n",mumps->id.ICNTL(30));
1679:       PetscViewerASCIIPrintf(viewer,"  ICNTL(31) (factors is discarded in the solve phase):    %d \n",mumps->id.ICNTL(31));
1680:       PetscViewerASCIIPrintf(viewer,"  ICNTL(33) (compute determinant):                        %d \n",mumps->id.ICNTL(33));

1682:       PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",mumps->id.CNTL(1));
1683:       PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2));
1684:       PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",mumps->id.CNTL(3));
1685:       PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",mumps->id.CNTL(4));
1686:       PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",mumps->id.CNTL(5));

1688:       /* infomation local to each processor */
1689:       PetscViewerASCIIPrintf(viewer, "  RINFO(1) (local estimated flops for the elimination after analysis): \n");
1690:       PetscViewerASCIIPushSynchronized(viewer);
1691:       PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %g \n",mumps->myid,mumps->id.RINFO(1));
1692:       PetscViewerFlush(viewer);
1693:       PetscViewerASCIIPrintf(viewer, "  RINFO(2) (local estimated flops for the assembly after factorization): \n");
1694:       PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(2));
1695:       PetscViewerFlush(viewer);
1696:       PetscViewerASCIIPrintf(viewer, "  RINFO(3) (local estimated flops for the elimination after factorization): \n");
1697:       PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(3));
1698:       PetscViewerFlush(viewer);

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

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

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

1712:       if (mumps->ninfo && mumps->ninfo <= 40){
1713:         PetscInt i;
1714:         for (i=0; i<mumps->ninfo; i++){
1715:           PetscViewerASCIIPrintf(viewer, "  INFO(%d): \n",mumps->info[i]);
1716:           PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(mumps->info[i]));
1717:           PetscViewerFlush(viewer);
1718:         }
1719:       }


1722:       PetscViewerASCIIPopSynchronized(viewer);

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

1730:         PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(3));
1731:         PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(4));
1732:         PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",mumps->id.INFOG(5));
1733:         PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",mumps->id.INFOG(6));
1734:         PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively use after analysis): %d \n",mumps->id.INFOG(7));
1735:         PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",mumps->id.INFOG(8));
1736:         PetscViewerASCIIPrintf(viewer,"  INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",mumps->id.INFOG(9));
1737:         PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",mumps->id.INFOG(10));
1738:         PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",mumps->id.INFOG(11));
1739:         PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",mumps->id.INFOG(12));
1740:         PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",mumps->id.INFOG(13));
1741:         PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",mumps->id.INFOG(14));
1742:         PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",mumps->id.INFOG(15));
1743:         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));
1744:         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));
1745:         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));
1746:         PetscViewerASCIIPrintf(viewer,"  INFOG(19) (size of all MUMPS internal data allocated during factorization: sum over all processors): %d \n",mumps->id.INFOG(19));
1747:         PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",mumps->id.INFOG(20));
1748:         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));
1749:         PetscViewerASCIIPrintf(viewer,"  INFOG(22) (size in MB of memory effectively used during factorization - sum over all processors): %d \n",mumps->id.INFOG(22));
1750:         PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",mumps->id.INFOG(23));
1751:         PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",mumps->id.INFOG(24));
1752:         PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",mumps->id.INFOG(25));
1753:         PetscViewerASCIIPrintf(viewer,"  INFOG(28) (after factorization: number of null pivots encountered): %d\n",mumps->id.INFOG(28));
1754:         PetscViewerASCIIPrintf(viewer,"  INFOG(29) (after factorization: effective number of entries in the factors (sum over all processors)): %d\n",mumps->id.INFOG(29));
1755:         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));
1756:         PetscViewerASCIIPrintf(viewer,"  INFOG(32) (after analysis: type of analysis done): %d\n",mumps->id.INFOG(32));
1757:         PetscViewerASCIIPrintf(viewer,"  INFOG(33) (value used for ICNTL(8)): %d\n",mumps->id.INFOG(33));
1758:         PetscViewerASCIIPrintf(viewer,"  INFOG(34) (exponent of the determinant if determinant is requested): %d\n",mumps->id.INFOG(34));
1759:       }
1760:     }
1761:   }
1762:   return(0);
1763: }

1765: PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
1766: {
1767:   Mat_MUMPS *mumps =(Mat_MUMPS*)A->data;

1770:   info->block_size        = 1.0;
1771:   info->nz_allocated      = mumps->id.INFOG(20);
1772:   info->nz_used           = mumps->id.INFOG(20);
1773:   info->nz_unneeded       = 0.0;
1774:   info->assemblies        = 0.0;
1775:   info->mallocs           = 0.0;
1776:   info->memory            = 0.0;
1777:   info->fill_ratio_given  = 0;
1778:   info->fill_ratio_needed = 0;
1779:   info->factor_mallocs    = 0;
1780:   return(0);
1781: }

1783: /* -------------------------------------------------------------------------------------------*/
1784: PetscErrorCode MatFactorSetSchurIS_MUMPS(Mat F, IS is)
1785: {
1786:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1787:   const PetscInt *idxs;
1788:   PetscInt       size,i;

1792:   if (mumps->size > 1) {
1793:     PetscBool ls,gs;

1795:     ISGetLocalSize(is,&size);
1796:     ls   = mumps->myid ? (size ? PETSC_FALSE : PETSC_TRUE) : PETSC_TRUE;
1797:     MPI_Allreduce(&ls,&gs,1,MPIU_BOOL,MPI_LAND,mumps->comm_mumps);
1798:     if (!gs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MUMPS distributed parallel Schur complements not yet supported from PETSc\n");
1799:   }
1800:   ISGetLocalSize(is,&size);
1801:   if (mumps->id.size_schur != size) {
1802:     PetscFree2(mumps->id.listvar_schur,mumps->id.schur);
1803:     mumps->id.size_schur = size;
1804:     mumps->id.schur_lld = size;
1805:     PetscMalloc2(size,&mumps->id.listvar_schur,size*size,&mumps->id.schur);
1806:   }
1807:   ISGetIndices(is,&idxs);
1808:   PetscMemcpy(mumps->id.listvar_schur,idxs,size*sizeof(PetscInt));
1809:   /* MUMPS expects Fortran style indices */
1810:   for (i=0;i<size;i++) mumps->id.listvar_schur[i]++;
1811:   ISRestoreIndices(is,&idxs);
1812:   if (mumps->size > 1) {
1813:     mumps->id.ICNTL(19) = 1; /* MUMPS returns Schur centralized on the host */
1814:   } else {
1815:     if (F->factortype == MAT_FACTOR_LU) {
1816:       mumps->id.ICNTL(19) = 3; /* MUMPS returns full matrix */
1817:     } else {
1818:       mumps->id.ICNTL(19) = 2; /* MUMPS returns lower triangular part */
1819:     }
1820:   }
1821:   /* set a special value of ICNTL (not handled my MUMPS) to be used in the solve phase by PETSc */
1822:   mumps->id.ICNTL(26) = -1;
1823:   return(0);
1824: }

1826: /* -------------------------------------------------------------------------------------------*/
1827: PetscErrorCode MatFactorCreateSchurComplement_MUMPS(Mat F,Mat* S)
1828: {
1829:   Mat            St;
1830:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1831:   PetscScalar    *array;
1832: #if defined(PETSC_USE_COMPLEX)
1833:   PetscScalar    im = PetscSqrtScalar((PetscScalar)-1.0);
1834: #endif

1838:   if (!mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it");
1839:   MatCreate(PETSC_COMM_SELF,&St);
1840:   MatSetSizes(St,PETSC_DECIDE,PETSC_DECIDE,mumps->id.size_schur,mumps->id.size_schur);
1841:   MatSetType(St,MATDENSE);
1842:   MatSetUp(St);
1843:   MatDenseGetArray(St,&array);
1844:   if (!mumps->sym) { /* MUMPS always return a full matrix */
1845:     if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
1846:       PetscInt i,j,N=mumps->id.size_schur;
1847:       for (i=0;i<N;i++) {
1848:         for (j=0;j<N;j++) {
1849: #if !defined(PETSC_USE_COMPLEX)
1850:           PetscScalar val = mumps->id.schur[i*N+j];
1851: #else
1852:           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
1853: #endif
1854:           array[j*N+i] = val;
1855:         }
1856:       }
1857:     } else { /* stored by columns */
1858:       PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));
1859:     }
1860:   } else { /* either full or lower-triangular (not packed) */
1861:     if (mumps->id.ICNTL(19) == 2) { /* lower triangular stored by columns */
1862:       PetscInt i,j,N=mumps->id.size_schur;
1863:       for (i=0;i<N;i++) {
1864:         for (j=i;j<N;j++) {
1865: #if !defined(PETSC_USE_COMPLEX)
1866:           PetscScalar val = mumps->id.schur[i*N+j];
1867: #else
1868:           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
1869: #endif
1870:           array[i*N+j] = val;
1871:           array[j*N+i] = val;
1872:         }
1873:       }
1874:     } else if (mumps->id.ICNTL(19) == 3) { /* full matrix */
1875:       PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));
1876:     } else { /* ICNTL(19) == 1 lower triangular stored by rows */
1877:       PetscInt i,j,N=mumps->id.size_schur;
1878:       for (i=0;i<N;i++) {
1879:         for (j=0;j<i+1;j++) {
1880: #if !defined(PETSC_USE_COMPLEX)
1881:           PetscScalar val = mumps->id.schur[i*N+j];
1882: #else
1883:           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
1884: #endif
1885:           array[i*N+j] = val;
1886:           array[j*N+i] = val;
1887:         }
1888:       }
1889:     }
1890:   }
1891:   MatDenseRestoreArray(St,&array);
1892:   *S = St;
1893:   return(0);
1894: }

1896: /* -------------------------------------------------------------------------------------------*/
1897: PetscErrorCode MatFactorGetSchurComplement_MUMPS(Mat F,Mat* S)
1898: {
1899:   Mat            St;
1900:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;

1904:   if (!mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it");
1905:   /* It should be the responsibility of the user to handle different ICNTL(19) cases and factorization stages if they want to work with the raw data */
1906:   MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.size_schur,(PetscScalar*)mumps->id.schur,&St);
1907:   *S = St;
1908:   return(0);
1909: }

1911: /* -------------------------------------------------------------------------------------------*/
1912: PetscErrorCode MatFactorFactorizeSchurComplement_MUMPS(Mat F)
1913: {
1914:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;

1918:   if (!mumps->id.ICNTL(19)) { /* do nothing */
1919:     return(0);
1920:   }
1921:   if (mumps->size > 1) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_SUP,"Parallel Schur complement not yet supported from PETSc");
1922:   MatMumpsFactorSchur_Private(mumps);
1923:   return(0);
1924: }

1926: PetscErrorCode MatFactorInvertSchurComplement_MUMPS(Mat F)
1927: {
1928:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;

1932:   if (!mumps->id.ICNTL(19)) { /* do nothing */
1933:     return(0);
1934:   }
1935:   if (mumps->size > 1) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_SUP,"Parallel Schur complement not yet supported from PETSc");
1936:   MatMumpsInvertSchur_Private(mumps);
1937:   return(0);
1938: }

1940: /* -------------------------------------------------------------------------------------------*/
1941: PetscErrorCode MatFactorSolveSchurComplement_MUMPS(Mat F, Vec rhs, Vec sol)
1942: {
1943:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1944:   MumpsScalar    *orhs;
1945:   PetscScalar    *osol,*nrhs,*nsol;
1946:   PetscInt       orhs_size,osol_size,olrhs_size;

1950:   if (!mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it");
1951:   if (mumps->size > 1) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_SUP,"Parallel Schur complement not yet supported from PETSc");

1953:   /* swap pointers */
1954:   orhs = mumps->id.redrhs;
1955:   olrhs_size = mumps->id.lredrhs;
1956:   orhs_size = mumps->sizeredrhs;
1957:   osol = mumps->schur_sol;
1958:   osol_size = mumps->schur_sizesol;
1959:   VecGetArray(rhs,&nrhs);
1960:   VecGetArray(sol,&nsol);
1961:   mumps->id.redrhs = (MumpsScalar*)nrhs;
1962:   VecGetLocalSize(rhs,&mumps->sizeredrhs);
1963:   mumps->id.lredrhs = mumps->sizeredrhs;
1964:   mumps->schur_sol = nsol;
1965:   VecGetLocalSize(sol,&mumps->schur_sizesol);

1967:   /* solve Schur complement */
1968:   mumps->id.nrhs = 1;
1969:   MatMumpsSolveSchur_Private(mumps,PETSC_FALSE);
1970:   /* restore pointers */
1971:   VecRestoreArray(rhs,&nrhs);
1972:   VecRestoreArray(sol,&nsol);
1973:   mumps->id.redrhs = orhs;
1974:   mumps->id.lredrhs = olrhs_size;
1975:   mumps->sizeredrhs = orhs_size;
1976:   mumps->schur_sol = osol;
1977:   mumps->schur_sizesol = osol_size;
1978:   return(0);
1979: }

1981: /* -------------------------------------------------------------------------------------------*/
1982: PetscErrorCode MatFactorSolveSchurComplementTranspose_MUMPS(Mat F, Vec rhs, Vec sol)
1983: {
1984:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1985:   MumpsScalar    *orhs;
1986:   PetscScalar    *osol,*nrhs,*nsol;
1987:   PetscInt       orhs_size,osol_size;

1991:   if (!mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it");
1992:   if (mumps->size > 1) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_SUP,"Parallel Schur complement not yet supported from PETSc");

1994:   /* swap pointers */
1995:   orhs = mumps->id.redrhs;
1996:   orhs_size = mumps->sizeredrhs;
1997:   osol = mumps->schur_sol;
1998:   osol_size = mumps->schur_sizesol;
1999:   VecGetArray(rhs,&nrhs);
2000:   VecGetArray(sol,&nsol);
2001:   mumps->id.redrhs = (MumpsScalar*)nrhs;
2002:   VecGetLocalSize(rhs,&mumps->sizeredrhs);
2003:   mumps->schur_sol = nsol;
2004:   VecGetLocalSize(sol,&mumps->schur_sizesol);

2006:   /* solve Schur complement */
2007:   mumps->id.nrhs = 1;
2008:   mumps->id.ICNTL(9) = 0;
2009:   MatMumpsSolveSchur_Private(mumps,PETSC_FALSE);
2010:   mumps->id.ICNTL(9) = 1;
2011:   /* restore pointers */
2012:   VecRestoreArray(rhs,&nrhs);
2013:   VecRestoreArray(sol,&nsol);
2014:   mumps->id.redrhs = orhs;
2015:   mumps->sizeredrhs = orhs_size;
2016:   mumps->schur_sol = osol;
2017:   mumps->schur_sizesol = osol_size;
2018:   return(0);
2019: }

2021: /* -------------------------------------------------------------------------------------------*/
2022: PetscErrorCode MatFactorSetSchurComplementSolverType_MUMPS(Mat F, PetscInt sym)
2023: {
2024:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

2027:   if (mumps->schur_factored && mumps->sym != mumps->schur_sym) {
2028:     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONG,"Cannot change the Schur solver! Schur complement data has been already factored");
2029:   }
2030:   mumps->schur_sym = sym;
2031:   return(0);
2032: }

2034: /* -------------------------------------------------------------------------------------------*/
2035: PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival)
2036: {
2037:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

2040:   mumps->id.ICNTL(icntl) = ival;
2041:   return(0);
2042: }

2044: PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt *ival)
2045: {
2046:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

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

2053: /*@
2054:   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()

2056:    Logically Collective on Mat

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

2063:   Options Database:
2064: .   -mat_mumps_icntl_<icntl> <ival>

2066:    Level: beginner

2068:    References:
2069: .     MUMPS Users' Guide

2071: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2072:  @*/
2073: PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival)
2074: {
2076: 
2079:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2082:   PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));
2083:   return(0);
2084: }

2086: /*@
2087:   MatMumpsGetIcntl - Get MUMPS parameter ICNTL()

2089:    Logically Collective on Mat

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

2095:   Output Parameter:
2096: .  ival - value of MUMPS ICNTL(icntl)

2098:    Level: beginner

2100:    References:
2101: .     MUMPS Users' Guide

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

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

2118: /* -------------------------------------------------------------------------------------------*/
2119: PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal val)
2120: {
2121:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

2124:   mumps->id.CNTL(icntl) = val;
2125:   return(0);
2126: }

2128: PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal *val)
2129: {
2130:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

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

2137: /*@
2138:   MatMumpsSetCntl - Set MUMPS parameter CNTL()

2140:    Logically Collective on Mat

2142:    Input Parameters:
2143: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2144: .  icntl - index of MUMPS parameter array CNTL()
2145: -  val - value of MUMPS CNTL(icntl)

2147:   Options Database:
2148: .   -mat_mumps_cntl_<icntl> <val>

2150:    Level: beginner

2152:    References:
2153: .     MUMPS Users' Guide

2155: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2156: @*/
2157: PetscErrorCode MatMumpsSetCntl(Mat F,PetscInt icntl,PetscReal val)
2158: {

2163:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2166:   PetscTryMethod(F,"MatMumpsSetCntl_C",(Mat,PetscInt,PetscReal),(F,icntl,val));
2167:   return(0);
2168: }

2170: /*@
2171:   MatMumpsGetCntl - Get MUMPS parameter CNTL()

2173:    Logically Collective on Mat

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

2179:   Output Parameter:
2180: .  val - value of MUMPS CNTL(icntl)

2182:    Level: beginner

2184:    References:
2185: .      MUMPS Users' Guide

2187: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2188: @*/
2189: PetscErrorCode MatMumpsGetCntl(Mat F,PetscInt icntl,PetscReal *val)
2190: {

2195:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2198:   PetscUseMethod(F,"MatMumpsGetCntl_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));
2199:   return(0);
2200: }

2202: PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F,PetscInt icntl,PetscInt *info)
2203: {
2204:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

2207:   *info = mumps->id.INFO(icntl);
2208:   return(0);
2209: }

2211: PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F,PetscInt icntl,PetscInt *infog)
2212: {
2213:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

2216:   *infog = mumps->id.INFOG(icntl);
2217:   return(0);
2218: }

2220: PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfo)
2221: {
2222:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

2225:   *rinfo = mumps->id.RINFO(icntl);
2226:   return(0);
2227: }

2229: PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfog)
2230: {
2231:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;

2234:   *rinfog = mumps->id.RINFOG(icntl);
2235:   return(0);
2236: }

2238: /*@
2239:   MatMumpsGetInfo - Get MUMPS parameter INFO()

2241:    Logically Collective on Mat

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

2247:   Output Parameter:
2248: .  ival - value of MUMPS INFO(icntl)

2250:    Level: beginner

2252:    References:
2253: .      MUMPS Users' Guide

2255: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2256: @*/
2257: PetscErrorCode MatMumpsGetInfo(Mat F,PetscInt icntl,PetscInt *ival)
2258: {

2263:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2265:   PetscUseMethod(F,"MatMumpsGetInfo_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));
2266:   return(0);
2267: }

2269: /*@
2270:   MatMumpsGetInfog - Get MUMPS parameter INFOG()

2272:    Logically Collective on Mat

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

2278:   Output Parameter:
2279: .  ival - value of MUMPS INFOG(icntl)

2281:    Level: beginner

2283:    References:
2284: .      MUMPS Users' Guide

2286: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2287: @*/
2288: PetscErrorCode MatMumpsGetInfog(Mat F,PetscInt icntl,PetscInt *ival)
2289: {

2294:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2296:   PetscUseMethod(F,"MatMumpsGetInfog_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));
2297:   return(0);
2298: }

2300: /*@
2301:   MatMumpsGetRinfo - Get MUMPS parameter RINFO()

2303:    Logically Collective on Mat

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

2309:   Output Parameter:
2310: .  val - value of MUMPS RINFO(icntl)

2312:    Level: beginner

2314:    References:
2315: .       MUMPS Users' Guide

2317: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2318: @*/
2319: PetscErrorCode MatMumpsGetRinfo(Mat F,PetscInt icntl,PetscReal *val)
2320: {

2325:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2327:   PetscUseMethod(F,"MatMumpsGetRinfo_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));
2328:   return(0);
2329: }

2331: /*@
2332:   MatMumpsGetRinfog - Get MUMPS parameter RINFOG()

2334:    Logically Collective on Mat

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

2340:   Output Parameter:
2341: .  val - value of MUMPS RINFOG(icntl)

2343:    Level: beginner

2345:    References:
2346: .      MUMPS Users' Guide

2348: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2349: @*/
2350: PetscErrorCode MatMumpsGetRinfog(Mat F,PetscInt icntl,PetscReal *val)
2351: {

2356:   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2358:   PetscUseMethod(F,"MatMumpsGetRinfog_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));
2359:   return(0);
2360: }

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

2366:   Works with MATAIJ and MATSBAIJ matrices

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

2370:   Use -pc_type cholesky or lu -pc_factor_mat_solver_package mumps to us this direct solver

2372:   Options Database Keys:
2373: +  -mat_mumps_icntl_1 - ICNTL(1): output stream for error messages 
2374: .  -mat_mumps_icntl_2 - ICNTL(2): output stream for diagnostic printing, statistics, and warning 
2375: .  -mat_mumps_icntl_3 -  ICNTL(3): output stream for global information, collected on the host 
2376: .  -mat_mumps_icntl_4 -  ICNTL(4): level of printing (0 to 4) 
2377: .  -mat_mumps_icntl_6 - ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7) 
2378: .  -mat_mumps_icntl_7 - ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis 
2379: .  -mat_mumps_icntl_8  - ICNTL(8): scaling strategy (-2 to 8 or 77) 
2380: .  -mat_mumps_icntl_10  - ICNTL(10): max num of refinements 
2381: .  -mat_mumps_icntl_11  - ICNTL(11): statistics related to an error analysis (via -ksp_view) 
2382: .  -mat_mumps_icntl_12  - ICNTL(12): an ordering strategy for symmetric matrices (0 to 3) 
2383: .  -mat_mumps_icntl_13  - ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting 
2384: .  -mat_mumps_icntl_14  - ICNTL(14): percentage increase in the estimated working space 
2385: .  -mat_mumps_icntl_19  - ICNTL(19): computes the Schur complement 
2386: .  -mat_mumps_icntl_22  - ICNTL(22): in-core/out-of-core factorization and solve (0 or 1) 
2387: .  -mat_mumps_icntl_23  - ICNTL(23): max size of the working memory (MB) that can allocate per processor 
2388: .  -mat_mumps_icntl_24  - ICNTL(24): detection of null pivot rows (0 or 1) 
2389: .  -mat_mumps_icntl_25  - ICNTL(25): compute a solution of a deficient matrix and a null space basis 
2390: .  -mat_mumps_icntl_26  - ICNTL(26): drives the solution phase if a Schur complement matrix 
2391: .  -mat_mumps_icntl_28  - ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering 
2392: .  -mat_mumps_icntl_29 - ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis 
2393: .  -mat_mumps_icntl_30 - ICNTL(30): compute user-specified set of entries in inv(A) 
2394: .  -mat_mumps_icntl_31 - ICNTL(31): indicates which factors may be discarded during factorization 
2395: .  -mat_mumps_icntl_33 - ICNTL(33): compute determinant 
2396: .  -mat_mumps_cntl_1  - CNTL(1): relative pivoting threshold 
2397: .  -mat_mumps_cntl_2  -  CNTL(2): stopping criterion of refinement 
2398: .  -mat_mumps_cntl_3 - CNTL(3): absolute pivoting threshold 
2399: .  -mat_mumps_cntl_4 - CNTL(4): value for static pivoting 
2400: -  -mat_mumps_cntl_5 - CNTL(5): fixation for null pivots 

2402:   Level: beginner

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

2411: .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage, MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog(), KSPGetPC(), PCGetFactor(), PCFactorGetMatrix()

2413: M*/

2415: static PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type)
2416: {
2418:   *type = MATSOLVERMUMPS;
2419:   return(0);
2420: }

2422: /* MatGetFactor for Seq and MPI AIJ matrices */
2423: static PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
2424: {
2425:   Mat            B;
2427:   Mat_MUMPS      *mumps;
2428:   PetscBool      isSeqAIJ;

2431:   /* Create the factorization matrix */
2432:   PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);
2433:   MatCreate(PetscObjectComm((PetscObject)A),&B);
2434:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
2435:   PetscStrallocpy("mumps",&((PetscObject)B)->type_name);
2436:   MatSetUp(B);

2438:   PetscNewLog(B,&mumps);

2440:   B->ops->view        = MatView_MUMPS;
2441:   B->ops->getinfo     = MatGetInfo_MUMPS;

2443:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);
2444:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);
2445:   PetscObjectComposeFunction((PetscObject)B,"MatFactorInvertSchurComplement_C",MatFactorInvertSchurComplement_MUMPS);
2446:   PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);
2447:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSchurComplement_C",MatFactorGetSchurComplement_MUMPS);
2448:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplement_C",MatFactorSolveSchurComplement_MUMPS);
2449:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplementTranspose_C",MatFactorSolveSchurComplementTranspose_MUMPS);
2450:   PetscObjectComposeFunction((PetscObject)B,"MatFactorFactorizeSchurComplement_C",MatFactorFactorizeSchurComplement_MUMPS);
2451:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurComplementSolverType_C",MatFactorSetSchurComplementSolverType_MUMPS);
2452:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);
2453:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);
2454:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);
2455:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);
2456:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);
2457:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);
2458:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);
2459:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);

2461:   if (ftype == MAT_FACTOR_LU) {
2462:     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
2463:     B->factortype            = MAT_FACTOR_LU;
2464:     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
2465:     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
2466:     mumps->sym = 0;
2467:   } else {
2468:     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2469:     B->factortype                  = MAT_FACTOR_CHOLESKY;
2470:     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
2471:     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
2472: #if defined(PETSC_USE_COMPLEX)
2473:     mumps->sym = 2;
2474: #else
2475:     if (A->spd_set && A->spd) mumps->sym = 1;
2476:     else                      mumps->sym = 2;
2477: #endif
2478:   }

2480:   /* set solvertype */
2481:   PetscFree(B->solvertype);
2482:   PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);

2484:   mumps->isAIJ    = PETSC_TRUE;
2485:   B->ops->destroy = MatDestroy_MUMPS;
2486:   B->data        = (void*)mumps;

2488:   PetscInitializeMUMPS(A,mumps);

2490:   *F = B;
2491:   return(0);
2492: }

2494: /* MatGetFactor for Seq and MPI SBAIJ matrices */
2495: static PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
2496: {
2497:   Mat            B;
2499:   Mat_MUMPS      *mumps;
2500:   PetscBool      isSeqSBAIJ;

2503:   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
2504:   if (A->rmap->bs > 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with block size > 1 with MUMPS Cholesky, use AIJ matrix instead");
2505:   PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);
2506:   /* Create the factorization matrix */
2507:   MatCreate(PetscObjectComm((PetscObject)A),&B);
2508:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
2509:   PetscStrallocpy("mumps",&((PetscObject)B)->type_name);
2510:   MatSetUp(B);

2512:   PetscNewLog(B,&mumps);
2513:   if (isSeqSBAIJ) {
2514:     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
2515:   } else {
2516:     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
2517:   }

2519:   B->ops->getinfo                = MatGetInfo_External;
2520:   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2521:   B->ops->view                   = MatView_MUMPS;

2523:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);
2524:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);
2525:   PetscObjectComposeFunction((PetscObject)B,"MatFactorInvertSchurComplement_C",MatFactorInvertSchurComplement_MUMPS);
2526:   PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);
2527:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSchurComplement_C",MatFactorGetSchurComplement_MUMPS);
2528:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplement_C",MatFactorSolveSchurComplement_MUMPS);
2529:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplementTranspose_C",MatFactorSolveSchurComplementTranspose_MUMPS);
2530:   PetscObjectComposeFunction((PetscObject)B,"MatFactorFactorizeSchurComplement_C",MatFactorFactorizeSchurComplement_MUMPS);
2531:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurComplementSolverType_C",MatFactorSetSchurComplementSolverType_MUMPS);
2532:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);
2533:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);
2534:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);
2535:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);
2536:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);
2537:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);
2538:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);
2539:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);

2541:   B->factortype = MAT_FACTOR_CHOLESKY;
2542: #if defined(PETSC_USE_COMPLEX)
2543:   mumps->sym = 2;
2544: #else
2545:   if (A->spd_set && A->spd) mumps->sym = 1;
2546:   else                      mumps->sym = 2;
2547: #endif

2549:   /* set solvertype */
2550:   PetscFree(B->solvertype);
2551:   PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);

2553:   mumps->isAIJ    = PETSC_FALSE;
2554:   B->ops->destroy = MatDestroy_MUMPS;
2555:   B->data        = (void*)mumps;

2557:   PetscInitializeMUMPS(A,mumps);

2559:   *F = B;
2560:   return(0);
2561: }

2563: static PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
2564: {
2565:   Mat            B;
2567:   Mat_MUMPS      *mumps;
2568:   PetscBool      isSeqBAIJ;

2571:   /* Create the factorization matrix */
2572:   PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);
2573:   MatCreate(PetscObjectComm((PetscObject)A),&B);
2574:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
2575:   PetscStrallocpy("mumps",&((PetscObject)B)->type_name);
2576:   MatSetUp(B);

2578:   PetscNewLog(B,&mumps);
2579:   if (ftype == MAT_FACTOR_LU) {
2580:     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
2581:     B->factortype            = MAT_FACTOR_LU;
2582:     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
2583:     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
2584:     mumps->sym = 0;
2585:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n");

2587:   B->ops->getinfo     = MatGetInfo_External;
2588:   B->ops->view        = MatView_MUMPS;

2590:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);
2591:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);
2592:   PetscObjectComposeFunction((PetscObject)B,"MatFactorInvertSchurComplement_C",MatFactorInvertSchurComplement_MUMPS);
2593:   PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);
2594:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSchurComplement_C",MatFactorGetSchurComplement_MUMPS);
2595:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplement_C",MatFactorSolveSchurComplement_MUMPS);
2596:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplementTranspose_C",MatFactorSolveSchurComplementTranspose_MUMPS);
2597:   PetscObjectComposeFunction((PetscObject)B,"MatFactorFactorizeSchurComplement_C",MatFactorFactorizeSchurComplement_MUMPS);
2598:   PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurComplementSolverType_C",MatFactorSetSchurComplementSolverType_MUMPS);
2599:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);
2600:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);
2601:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);
2602:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);
2603:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);
2604:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);
2605:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);
2606:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);

2608:   /* set solvertype */
2609:   PetscFree(B->solvertype);
2610:   PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);

2612:   mumps->isAIJ    = PETSC_TRUE;
2613:   B->ops->destroy = MatDestroy_MUMPS;
2614:   B->data        = (void*)mumps;

2616:   PetscInitializeMUMPS(A,mumps);

2618:   *F = B;
2619:   return(0);
2620: }

2622: PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_MUMPS(void)
2623: {

2627:   MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);
2628:   MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);
2629:   MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);
2630:   MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);
2631:   MatSolverPackageRegister(MATSOLVERMUMPS,MATMPISBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);
2632:   MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);
2633:   MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);
2634:   MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);
2635:   MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);
2636:   MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);
2637:   return(0);
2638: }