Actual source code: mumps.c

petsc-3.6.1 2015-08-06
Report Typos and Errors
  2: /*
  3:     Provides an interface to the MUMPS sparse solver
  4: */

  6: #include <../src/mat/impls/aij/mpi/mpiaij.h> /*I  "petscmat.h"  I*/
  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,CleanUpMUMPS;
 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:   PetscBool    schur_second_solve;
 91:   PetscInt     sizeredrhs;
 92:   PetscInt     *schur_pivots;
 93:   PetscInt     schur_B_lwork;
 94:   PetscScalar  *schur_work;
 95:   PetscScalar  *schur_sol;
 96:   PetscInt     schur_sizesol;
 97:   PetscBool    schur_restored;
 98:   PetscBool    schur_factored;
 99:   PetscBool    schur_inverted;

101:   PetscErrorCode (*Destroy)(Mat);
102:   PetscErrorCode (*ConvertToTriples)(Mat, int, MatReuse, int*, int**, int**, PetscScalar**);
103: } Mat_MUMPS;

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

109: static PetscErrorCode MatMumpsResetSchur_Private(Mat_MUMPS* mumps)
110: {

114:   PetscFree2(mumps->id.listvar_schur,mumps->id.schur);
115:   PetscFree(mumps->id.redrhs);
116:   PetscFree(mumps->schur_sol);
117:   PetscFree(mumps->schur_pivots);
118:   PetscFree(mumps->schur_work);
119:   if (!mumps->schur_restored) {
120:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Schur matrix has not been restored using MatMumpsRestoreSchurComplement");
121:   }
122:   mumps->id.size_schur = 0;
123:   mumps->id.ICNTL(19) = 0;
124:   return(0);
125: }

129: static PetscErrorCode MatMumpsFactorSchur_Private(Mat_MUMPS* mumps)
130: {
131:   PetscBLASInt   B_N,B_ierr,B_slda;

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

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

185: static PetscErrorCode MatMumpsInvertSchur_Private(Mat_MUMPS* mumps)
186: {
187:   PetscBLASInt   B_N,B_ierr,B_slda;

191:   MatMumpsFactorSchur_Private(mumps);
192:   PetscBLASIntCast(mumps->id.size_schur,&B_N);
193:   PetscBLASIntCast(mumps->id.schur_lld,&B_slda);
194:   if (!mumps->sym) { /* MUMPS always return a full Schur matrix */
195:     if (!mumps->schur_work) {
196:       PetscScalar lwork;

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

235: static PetscErrorCode MatMumpsSolveSchur_Private(Mat_MUMPS* mumps, PetscBool sol_in_redrhs)
236: {
237:   PetscBLASInt   B_N,B_Nrhs,B_ierr,B_slda,B_rlda;
238:   PetscScalar    one=1.,zero=0.;

242:   MatMumpsFactorSchur_Private(mumps);
243:   PetscBLASIntCast(mumps->id.size_schur,&B_N);
244:   PetscBLASIntCast(mumps->id.schur_lld,&B_slda);
245:   PetscBLASIntCast(mumps->id.nrhs,&B_Nrhs);
246:   PetscBLASIntCast(mumps->id.lredrhs,&B_rlda);
247:   if (mumps->schur_inverted) {
248:     PetscInt sizesol = B_Nrhs*B_N;
249:     if (!mumps->schur_sol || sizesol > mumps->schur_sizesol) {
250:       PetscFree(mumps->schur_sol);
251:       PetscMalloc1(sizesol,&mumps->schur_sol);
252:       mumps->schur_sizesol = sizesol;
253:     }
254:     if (!mumps->sym) {
255:       char type[2];
256:       if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
257:         if (!mumps->id.ICNTL(9)) { /* transpose solve */
258:           sprintf(type,"N");
259:         } else {
260:           sprintf(type,"T");
261:         }
262:       } else { /* stored by columns */
263:         if (!mumps->id.ICNTL(9)) { /* transpose solve */
264:           sprintf(type,"T");
265:         } else {
266:           sprintf(type,"N");
267:         }
268:       }
269:       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));
270:     } else {
271:       char ord[2];
272:       if (mumps->id.ICNTL(19) == 2 || mumps->id.ICNTL(19) == 3) { /* lower triangular stored by columns or full matrix */
273:         sprintf(ord,"L");
274:       } else { /* ICNTL(19) == 1 lower triangular stored by rows */
275:         sprintf(ord,"U");
276:       }
277:       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));
278:     }
279:     if (sol_in_redrhs) {
280:       PetscMemcpy(mumps->id.redrhs,mumps->schur_sol,sizesol*sizeof(PetscScalar));
281:     }
282:   } else { /* Schur complement has not been inverted */
283:     MumpsScalar *orhs=NULL;

285:     if (!sol_in_redrhs) {
286:       PetscInt sizesol = B_Nrhs*B_N;
287:       if (!mumps->schur_sol || sizesol > mumps->schur_sizesol) {
288:         PetscFree(mumps->schur_sol);
289:         PetscMalloc1(sizesol,&mumps->schur_sol);
290:         mumps->schur_sizesol = sizesol;
291:       }
292:       orhs = mumps->id.redrhs;
293:       PetscMemcpy(mumps->schur_sol,mumps->id.redrhs,sizesol*sizeof(PetscScalar));
294:       mumps->id.redrhs = (MumpsScalar*)mumps->schur_sol;
295:     }
296:     if (!mumps->sym) { /* MUMPS always return a full Schur matrix */
297:       char type[2];
298:       if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
299:         if (!mumps->id.ICNTL(9)) { /* transpose solve */
300:           sprintf(type,"N");
301:         } else {
302:           sprintf(type,"T");
303:         }
304:       } else { /* stored by columns */
305:         if (!mumps->id.ICNTL(9)) { /* transpose solve */
306:           sprintf(type,"T");
307:         } else {
308:           sprintf(type,"N");
309:         }
310:       }
311:       PetscFPTrapPush(PETSC_FP_TRAP_OFF);
312:       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));
313:       PetscFPTrapPop();
314:       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRS Lapack routine %d",(int)B_ierr);
315:     } else { /* either full or lower-triangular (not packed) */
316:       char ord[2];
317:       if (mumps->id.ICNTL(19) == 2 || mumps->id.ICNTL(19) == 3) { /* lower triangular stored by columns or full matrix */
318:         sprintf(ord,"L");
319:       } else { /* ICNTL(19) == 1 lower triangular stored by rows */
320:         sprintf(ord,"U");
321:       }
322:       if (mumps->id.sym == 2) {
323:         PetscFPTrapPush(PETSC_FP_TRAP_OFF);
324:         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));
325:         PetscFPTrapPop();
326:         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRS Lapack routine %d",(int)B_ierr);
327:       } else {
328:         PetscFPTrapPush(PETSC_FP_TRAP_OFF);
329:         PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_(ord,&B_N,&B_Nrhs,(PetscScalar*)mumps->id.schur,&B_slda,(PetscScalar*)mumps->id.redrhs,&B_rlda,&B_ierr));
330:         PetscFPTrapPop();
331:         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRS Lapack routine %d",(int)B_ierr);
332:       }
333:     }
334:     if (!sol_in_redrhs) {
335:       mumps->id.redrhs = orhs;
336:     }
337:   }
338:   return(0);
339: }

343: static PetscErrorCode MatMumpsHandleSchur_Private(Mat_MUMPS* mumps)
344: {

348:   if (!mumps->id.ICNTL(19)) { /* do nothing when Schur complement has not been computed */
349:     return(0);
350:   }
351:   if (!mumps->schur_second_solve) { /* prepare for the condensation step */
352:     /* check if schur complement has been computed
353:        We set by default ICNTL(26) == -1 when Schur indices have been provided by the user.
354:        According to MUMPS (5.0.0) manual, any value should be harmful during the factorization phase
355:        Unless the user provides a valid value for ICNTL(26), MatSolve and MatMatSolve routines solve the full system.
356:        This requires an extra call to PetscMUMPS_c and the computation of the factors for S, handled setting double_schur_solve to PETSC_TRUE */
357:     if (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2) {
358:       PetscInt sizeredrhs = mumps->id.nrhs*mumps->id.size_schur;
359:       /* allocate MUMPS internal array to store reduced right-hand sides */
360:       if (!mumps->id.redrhs || sizeredrhs > mumps->sizeredrhs) {
361:         PetscFree(mumps->id.redrhs);
362:         mumps->id.lredrhs = mumps->id.size_schur;
363:         PetscMalloc1(mumps->id.nrhs*mumps->id.lredrhs,&mumps->id.redrhs);
364:         mumps->sizeredrhs = mumps->id.nrhs*mumps->id.lredrhs;
365:       }
366:       mumps->schur_second_solve = PETSC_TRUE;
367:       mumps->id.ICNTL(26) = 1; /* condensation phase */
368:     }
369:   } else { /* prepare for the expansion step */
370:     /* solve Schur complement (this should be done by the MUMPS user, so basically us) */
371:     MatMumpsSolveSchur_Private(mumps,PETSC_TRUE);
372:     mumps->id.ICNTL(26) = 2; /* expansion phase */
373:     PetscMUMPS_c(&mumps->id);
374:     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));
375:     /* restore defaults */
376:     mumps->id.ICNTL(26) = -1;
377:     mumps->schur_second_solve = PETSC_FALSE;
378:   }
379:   return(0);
380: }

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

385:   input:
386:     A       - matrix in aij,baij or sbaij (bs=1) format
387:     shift   - 0: C style output triple; 1: Fortran style output triple.
388:     reuse   - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple
389:               MAT_REUSE_MATRIX:   only the values in v array are updated
390:   output:
391:     nnz     - dim of r, c, and v (number of local nonzero entries of A)
392:     r, c, v - row and col index, matrix values (matrix triples)

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

398:  */

402: PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
403: {
404:   const PetscInt *ai,*aj,*ajj,M=A->rmap->n;
405:   PetscInt       nz,rnz,i,j;
407:   PetscInt       *row,*col;
408:   Mat_SeqAIJ     *aa=(Mat_SeqAIJ*)A->data;

411:   *v=aa->a;
412:   if (reuse == MAT_INITIAL_MATRIX) {
413:     nz   = aa->nz;
414:     ai   = aa->i;
415:     aj   = aa->j;
416:     *nnz = nz;
417:     PetscMalloc1(2*nz, &row);
418:     col  = row + nz;

420:     nz = 0;
421:     for (i=0; i<M; i++) {
422:       rnz = ai[i+1] - ai[i];
423:       ajj = aj + ai[i];
424:       for (j=0; j<rnz; j++) {
425:         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
426:       }
427:     }
428:     *r = row; *c = col;
429:   }
430:   return(0);
431: }

435: PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
436: {
437:   Mat_SeqBAIJ    *aa=(Mat_SeqBAIJ*)A->data;
438:   const PetscInt *ai,*aj,*ajj,bs2 = aa->bs2;
439:   PetscInt       bs,M,nz,idx=0,rnz,i,j,k,m;
441:   PetscInt       *row,*col;

444:   MatGetBlockSize(A,&bs);
445:   M = A->rmap->N/bs;
446:   *v = aa->a;
447:   if (reuse == MAT_INITIAL_MATRIX) {
448:     ai   = aa->i; aj = aa->j;
449:     nz   = bs2*aa->nz;
450:     *nnz = nz;
451:     PetscMalloc1(2*nz, &row);
452:     col  = row + nz;

454:     for (i=0; i<M; i++) {
455:       ajj = aj + ai[i];
456:       rnz = ai[i+1] - ai[i];
457:       for (k=0; k<rnz; k++) {
458:         for (j=0; j<bs; j++) {
459:           for (m=0; m<bs; m++) {
460:             row[idx]   = i*bs + m + shift;
461:             col[idx++] = bs*(ajj[k]) + j + shift;
462:           }
463:         }
464:       }
465:     }
466:     *r = row; *c = col;
467:   }
468:   return(0);
469: }

473: PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
474: {
475:   const PetscInt *ai, *aj,*ajj,M=A->rmap->n;
476:   PetscInt       nz,rnz,i,j;
478:   PetscInt       *row,*col;
479:   Mat_SeqSBAIJ   *aa=(Mat_SeqSBAIJ*)A->data;

482:   *v = aa->a;
483:   if (reuse == MAT_INITIAL_MATRIX) {
484:     nz   = aa->nz;
485:     ai   = aa->i;
486:     aj   = aa->j;
487:     *v   = aa->a;
488:     *nnz = nz;
489:     PetscMalloc1(2*nz, &row);
490:     col  = row + nz;

492:     nz = 0;
493:     for (i=0; i<M; i++) {
494:       rnz = ai[i+1] - ai[i];
495:       ajj = aj + ai[i];
496:       for (j=0; j<rnz; j++) {
497:         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
498:       }
499:     }
500:     *r = row; *c = col;
501:   }
502:   return(0);
503: }

507: PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
508: {
509:   const PetscInt    *ai,*aj,*ajj,*adiag,M=A->rmap->n;
510:   PetscInt          nz,rnz,i,j;
511:   const PetscScalar *av,*v1;
512:   PetscScalar       *val;
513:   PetscErrorCode    ierr;
514:   PetscInt          *row,*col;
515:   Mat_SeqAIJ        *aa=(Mat_SeqAIJ*)A->data;

518:   ai   =aa->i; aj=aa->j;av=aa->a;
519:   adiag=aa->diag;
520:   if (reuse == MAT_INITIAL_MATRIX) {
521:     /* count nz in the uppper triangular part of A */
522:     nz = 0;
523:     for (i=0; i<M; i++) nz += ai[i+1] - adiag[i];
524:     *nnz = nz;

526:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
527:     col  = row + nz;
528:     val  = (PetscScalar*)(col + nz);

530:     nz = 0;
531:     for (i=0; i<M; i++) {
532:       rnz = ai[i+1] - adiag[i];
533:       ajj = aj + adiag[i];
534:       v1  = av + adiag[i];
535:       for (j=0; j<rnz; j++) {
536:         row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
537:       }
538:     }
539:     *r = row; *c = col; *v = val;
540:   } else {
541:     nz = 0; val = *v;
542:     for (i=0; i <M; i++) {
543:       rnz = ai[i+1] - adiag[i];
544:       ajj = aj + adiag[i];
545:       v1  = av + adiag[i];
546:       for (j=0; j<rnz; j++) {
547:         val[nz++] = v1[j];
548:       }
549:     }
550:   }
551:   return(0);
552: }

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

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

572:   garray = mat->garray;

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

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

586:   jj = 0; irow = rstart;
587:   for (i=0; i<m; i++) {
588:     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
589:     countA = ai[i+1] - ai[i];
590:     countB = bi[i+1] - bi[i];
591:     bjj    = bj + bi[i];
592:     v1     = av + ai[i];
593:     v2     = bv + bi[i];

595:     /* A-part */
596:     for (j=0; j<countA; j++) {
597:       if (reuse == MAT_INITIAL_MATRIX) {
598:         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
599:       }
600:       val[jj++] = v1[j];
601:     }

603:     /* B-part */
604:     for (j=0; j < countB; j++) {
605:       if (reuse == MAT_INITIAL_MATRIX) {
606:         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
607:       }
608:       val[jj++] = v2[j];
609:     }
610:     irow++;
611:   }
612:   return(0);
613: }

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

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

633:   garray = mat->garray;

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

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

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

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

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

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

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

701:     *r = row; *c = col; *v = val;
702:   } else {
703:     row = *r; col = *c; val = *v;
704:   }

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

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

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

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

762:   ai=aa->i; aj=aa->j; adiag=aa->diag;
763:   bi=bb->i; bj=bb->j; garray = mat->garray;
764:   av=aa->a; bv=bb->a;

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

768:   if (reuse == MAT_INITIAL_MATRIX) {
769:     nza = 0;    /* num of upper triangular entries in mat->A, including diagonals */
770:     nzb = 0;    /* num of upper triangular entries in mat->B */
771:     for (i=0; i<m; i++) {
772:       nza   += (ai[i+1] - adiag[i]);
773:       countB = bi[i+1] - bi[i];
774:       bjj    = bj + bi[i];
775:       for (j=0; j<countB; j++) {
776:         if (garray[bjj[j]] > rstart) nzb++;
777:       }
778:     }

780:     nz   = nza + nzb; /* total nz of upper triangular part of mat */
781:     *nnz = nz;
782:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
783:     col  = row + nz;
784:     val  = (PetscScalar*)(col + nz);

786:     *r = row; *c = col; *v = val;
787:   } else {
788:     row = *r; col = *c; val = *v;
789:   }

791:   jj = 0; irow = rstart;
792:   for (i=0; i<m; i++) {
793:     ajj    = aj + adiag[i];                 /* ptr to the beginning of the diagonal of this row */
794:     v1     = av + adiag[i];
795:     countA = ai[i+1] - adiag[i];
796:     countB = bi[i+1] - bi[i];
797:     bjj    = bj + bi[i];
798:     v2     = bv + bi[i];

800:     /* A-part */
801:     for (j=0; j<countA; j++) {
802:       if (reuse == MAT_INITIAL_MATRIX) {
803:         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
804:       }
805:       val[jj++] = v1[j];
806:     }

808:     /* B-part */
809:     for (j=0; j < countB; j++) {
810:       if (garray[bjj[j]] > rstart) {
811:         if (reuse == MAT_INITIAL_MATRIX) {
812:           row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
813:         }
814:         val[jj++] = v2[j];
815:       }
816:     }
817:     irow++;
818:   }
819:   return(0);
820: }

824: PetscErrorCode MatGetDiagonal_MUMPS(Mat A,Vec v)
825: {
827:   SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type: MUMPS factor");
828:   return(0);
829: }

833: PetscErrorCode MatDestroy_MUMPS(Mat A)
834: {
835:   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;

839:   if (mumps->CleanUpMUMPS) {
840:     /* Terminate instance, deallocate memories */
841:     PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);
842:     VecScatterDestroy(&mumps->scat_rhs);
843:     VecScatterDestroy(&mumps->scat_sol);
844:     VecDestroy(&mumps->b_seq);
845:     VecDestroy(&mumps->x_seq);
846:     PetscFree(mumps->id.perm_in);
847:     PetscFree(mumps->irn);
848:     PetscFree(mumps->info);
849:     MatMumpsResetSchur_Private(mumps);
850:     mumps->id.job = JOB_END;
851:     PetscMUMPS_c(&mumps->id);
852:     MPI_Comm_free(&(mumps->comm_mumps));
853:   }
854:   if (mumps->Destroy) {
855:     (mumps->Destroy)(A);
856:   }
857:   PetscFree(A->spptr);

859:   /* clear composed functions */
860:   PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);
861:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetIcntl_C",NULL);
862:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetIcntl_C",NULL);
863:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetCntl_C",NULL);
864:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetCntl_C",NULL);

866:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfo_C",NULL);
867:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfog_C",NULL);
868:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfo_C",NULL);
869:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfog_C",NULL);

871:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetSchurIndices_C",NULL);
872:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsInvertSchurComplement_C",NULL);
873:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsCreateSchurComplement_C",NULL);
874:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetSchurComplement_C",NULL);
875:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsRestoreSchurComplement_C",NULL);
876:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsSolveSchurComplement_C",NULL);
877:   PetscObjectComposeFunction((PetscObject)A,"MatMumpsSolveSchurComplementTranspose_C",NULL);
878:   return(0);
879: }

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

894:   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);
895:   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:   mumps->id.nrhs = 1;
897:   b_seq          = mumps->b_seq;
898:   if (mumps->size > 1) {
899:     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
900:     VecScatterBegin(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);
901:     VecScatterEnd(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);
902:     if (!mumps->myid) {VecGetArray(b_seq,&array);}
903:   } else {  /* size == 1 */
904:     VecCopy(b,x);
905:     VecGetArray(x,&array);
906:   }
907:   if (!mumps->myid) { /* define rhs on the host */
908:     mumps->id.nrhs = 1;
909:     mumps->id.rhs = (MumpsScalar*)array;
910:   }

912:   /* handle condensation step of Schur complement (if any) */
913:   MatMumpsHandleSchur_Private(mumps);

915:   /* solve phase */
916:   /*-------------*/
917:   mumps->id.job = JOB_SOLVE;
918:   PetscMUMPS_c(&mumps->id);
919:   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));

921:   /* handle expansion step of Schur complement (if any) */
922:   MatMumpsHandleSchur_Private(mumps);

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

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

942:     VecScatterBegin(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);
943:     VecScatterEnd(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);
944:   }
945:   return(0);
946: }

950: PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x)
951: {
952:   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;

956:   mumps->id.ICNTL(9) = 0;
957:   MatSolve_MUMPS(A,b,x);
958:   mumps->id.ICNTL(9) = 1;
959:   return(0);
960: }

964: PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X)
965: {
967:   PetscBool      flg;
968:   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
969:   PetscInt       i,nrhs,M;
970:   PetscScalar    *array,*bray;

973:   PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
974:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
975:   PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
976:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
977:   if (B->rmap->n != X->rmap->n) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B and X must have same row distribution");

979:   MatGetSize(B,&M,&nrhs);
980:   mumps->id.nrhs = nrhs;
981:   mumps->id.lrhs = M;

983:   if (mumps->size == 1) {
984:     /* copy B to X */
985:     MatDenseGetArray(B,&bray);
986:     MatDenseGetArray(X,&array);
987:     PetscMemcpy(array,bray,M*nrhs*sizeof(PetscScalar));
988:     MatDenseRestoreArray(B,&bray);
989:     mumps->id.rhs = (MumpsScalar*)array;
990:     /* handle condensation step of Schur complement (if any) */
991:     MatMumpsHandleSchur_Private(mumps);

993:     /* solve phase */
994:     /*-------------*/
995:     mumps->id.job = JOB_SOLVE;
996:     PetscMUMPS_c(&mumps->id);
997:     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));

999:     /* handle expansion step of Schur complement (if any) */
1000:     MatMumpsHandleSchur_Private(mumps);
1001:     MatDenseRestoreArray(X,&array);
1002:   } else {  /*--------- parallel case --------*/
1003:     PetscInt       lsol_loc,nlsol_loc,*isol_loc,*idx,*iidx,*idxx,*isol_loc_save;
1004:     MumpsScalar    *sol_loc,*sol_loc_save;
1005:     IS             is_to,is_from;
1006:     PetscInt       k,proc,j,m;
1007:     const PetscInt *rstart;
1008:     Vec            v_mpi,b_seq,x_seq;
1009:     VecScatter     scat_rhs,scat_sol;

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

1015:     lsol_loc  = mumps->id.INFO(23);
1016:     nlsol_loc = nrhs*lsol_loc;     /* length of sol_loc */
1017:     PetscMalloc2(nlsol_loc,&sol_loc,nlsol_loc,&isol_loc);
1018:     mumps->id.sol_loc = (MumpsScalar*)sol_loc;
1019:     mumps->id.isol_loc = isol_loc;

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

1023:     /* copy rhs matrix B into vector v_mpi */
1024:     MatGetLocalSize(B,&m,NULL);
1025:     MatDenseGetArray(B,&bray);
1026:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)B),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi);
1027:     MatDenseRestoreArray(B,&bray);

1029:     /* scatter v_mpi to b_seq because MUMPS only supports centralized rhs */
1030:     /* idx: maps from k-th index of v_mpi to (i,j)-th global entry of B;
1031:       iidx: inverse of idx, will be used by scattering xx_seq -> X       */
1032:     PetscMalloc2(nrhs*M,&idx,nrhs*M,&iidx);
1033:     MatGetOwnershipRanges(B,&rstart);
1034:     k = 0;
1035:     for (proc=0; proc<mumps->size; proc++){
1036:       for (j=0; j<nrhs; j++){
1037:         for (i=rstart[proc]; i<rstart[proc+1]; i++){
1038:           iidx[j*M + i] = k;
1039:           idx[k++]      = j*M + i;
1040:         }
1041:       }
1042:     }

1044:     if (!mumps->myid) {
1045:       VecCreateSeq(PETSC_COMM_SELF,nrhs*M,&b_seq);
1046:       ISCreateGeneral(PETSC_COMM_SELF,nrhs*M,idx,PETSC_COPY_VALUES,&is_to);
1047:       ISCreateStride(PETSC_COMM_SELF,nrhs*M,0,1,&is_from);
1048:     } else {
1049:       VecCreateSeq(PETSC_COMM_SELF,0,&b_seq);
1050:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_to);
1051:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_from);
1052:     }
1053:     VecScatterCreate(v_mpi,is_from,b_seq,is_to,&scat_rhs);
1054:     VecScatterBegin(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);
1055:     ISDestroy(&is_to);
1056:     ISDestroy(&is_from);
1057:     VecScatterEnd(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);

1059:     if (!mumps->myid) { /* define rhs on the host */
1060:       VecGetArray(b_seq,&bray);
1061:       mumps->id.rhs = (MumpsScalar*)bray;
1062:       VecRestoreArray(b_seq,&bray);
1063:     }

1065:     /* solve phase */
1066:     /*-------------*/
1067:     mumps->id.job = JOB_SOLVE;
1068:     PetscMUMPS_c(&mumps->id);
1069:     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));

1071:     /* scatter mumps distributed solution to petsc vector v_mpi, which shares local arrays with solution matrix X */
1072:     MatDenseGetArray(X,&array);
1073:     VecPlaceArray(v_mpi,array);
1074: 
1075:     /* create scatter scat_sol */
1076:     PetscMalloc1(nlsol_loc,&idxx);
1077:     ISCreateStride(PETSC_COMM_SELF,nlsol_loc,0,1,&is_from);
1078:     for (i=0; i<lsol_loc; i++) {
1079:       isol_loc[i] -= 1; /* change Fortran style to C style */
1080:       idxx[i] = iidx[isol_loc[i]];
1081:       for (j=1; j<nrhs; j++){
1082:         idxx[j*lsol_loc+i] = iidx[isol_loc[i]+j*M];
1083:       }
1084:     }
1085:     ISCreateGeneral(PETSC_COMM_SELF,nlsol_loc,idxx,PETSC_COPY_VALUES,&is_to);
1086:     VecScatterCreate(x_seq,is_from,v_mpi,is_to,&scat_sol);
1087:     VecScatterBegin(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);
1088:     ISDestroy(&is_from);
1089:     ISDestroy(&is_to);
1090:     VecScatterEnd(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);
1091:     MatDenseRestoreArray(X,&array);

1093:     /* free spaces */
1094:     mumps->id.sol_loc = sol_loc_save;
1095:     mumps->id.isol_loc = isol_loc_save;

1097:     PetscFree2(sol_loc,isol_loc);
1098:     PetscFree2(idx,iidx);
1099:     PetscFree(idxx);
1100:     VecDestroy(&x_seq);
1101:     VecDestroy(&v_mpi);
1102:     VecDestroy(&b_seq);
1103:     VecScatterDestroy(&scat_rhs);
1104:     VecScatterDestroy(&scat_sol);
1105:   }
1106:   return(0);
1107: }

1109: #if !defined(PETSC_USE_COMPLEX)
1110: /*
1111:   input:
1112:    F:        numeric factor
1113:   output:
1114:    nneg:     total number of negative pivots
1115:    nzero:    0
1116:    npos:     (global dimension of F) - nneg
1117: */

1121: PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
1122: {
1123:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
1125:   PetscMPIInt    size;

1128:   MPI_Comm_size(PetscObjectComm((PetscObject)F),&size);
1129:   /* 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 */
1130:   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));

1132:   if (nneg) *nneg = mumps->id.INFOG(12);
1133:   if (nzero || npos) {
1134:     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");
1135:     if (nzero) *nzero = mumps->id.INFOG(28);
1136:     if (npos) *npos   = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28));
1137:   }
1138:   return(0);
1139: }
1140: #endif /* !defined(PETSC_USE_COMPLEX) */

1144: PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
1145: {
1146:   Mat_MUMPS      *mumps =(Mat_MUMPS*)(F)->spptr;
1148:   Mat            F_diag;
1149:   PetscBool      isMPIAIJ;

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

1154:   /* numerical factorization phase */
1155:   /*-------------------------------*/
1156:   mumps->id.job = JOB_FACTNUMERIC;
1157:   if (!mumps->id.ICNTL(18)) { /* A is centralized */
1158:     if (!mumps->myid) {
1159:       mumps->id.a = (MumpsScalar*)mumps->val;
1160:     }
1161:   } else {
1162:     mumps->id.a_loc = (MumpsScalar*)mumps->val;
1163:   }
1164:   PetscMUMPS_c(&mumps->id);
1165:   if (mumps->id.INFOG(1) < 0) {
1166:     if (mumps->id.INFO(1) == -13) {
1167:       if (mumps->id.INFO(2) < 0) {
1168:         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",-mumps->id.INFO(2));
1169:       } else {
1170:         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d bytes\n",mumps->id.INFO(2));
1171:       }
1172:     } else SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: INFO(1)=%d, INFO(2)=%d\n",mumps->id.INFO(1),mumps->id.INFO(2));
1173:   }
1174:   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));
1175: 
1176:   (F)->assembled        = PETSC_TRUE;
1177:   mumps->matstruc       = SAME_NONZERO_PATTERN;
1178:   mumps->CleanUpMUMPS   = PETSC_TRUE;
1179:   mumps->schur_factored = PETSC_FALSE;
1180:   mumps->schur_inverted = PETSC_FALSE;

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

1185:   if (mumps->size > 1) {
1186:     PetscInt    lsol_loc;
1187:     PetscScalar *sol_loc;

1189:     PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);
1190:     if (isMPIAIJ) F_diag = ((Mat_MPIAIJ*)(F)->data)->A;
1191:     else F_diag = ((Mat_MPISBAIJ*)(F)->data)->A;
1192:     F_diag->assembled = PETSC_TRUE;

1194:     /* distributed solution; Create x_seq=sol_loc for repeated use */
1195:     if (mumps->x_seq) {
1196:       VecScatterDestroy(&mumps->scat_sol);
1197:       PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);
1198:       VecDestroy(&mumps->x_seq);
1199:     }
1200:     lsol_loc = mumps->id.INFO(23); /* length of sol_loc */
1201:     PetscMalloc2(lsol_loc,&sol_loc,lsol_loc,&mumps->id.isol_loc);
1202:     mumps->id.lsol_loc = lsol_loc;
1203:     mumps->id.sol_loc = (MumpsScalar*)sol_loc;
1204:     VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq);
1205:   }
1206:   return(0);
1207: }

1209: /* Sets MUMPS options from the options database */
1212: PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A)
1213: {
1214:   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
1216:   PetscInt       icntl,info[40],i,ninfo=40;
1217:   PetscBool      flg;

1220:   PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MUMPS Options","Mat");
1221:   PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);
1222:   if (flg) mumps->id.ICNTL(1) = icntl;
1223:   PetscOptionsInt("-mat_mumps_icntl_2","ICNTL(2): output stream for diagnostic printing, statistics, and warning","None",mumps->id.ICNTL(2),&icntl,&flg);
1224:   if (flg) mumps->id.ICNTL(2) = icntl;
1225:   PetscOptionsInt("-mat_mumps_icntl_3","ICNTL(3): output stream for global information, collected on the host","None",mumps->id.ICNTL(3),&icntl,&flg);
1226:   if (flg) mumps->id.ICNTL(3) = icntl;

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

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

1235:   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);
1236:   if (flg) {
1237:     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");
1238:     else mumps->id.ICNTL(7) = icntl;
1239:   }

1241:   PetscOptionsInt("-mat_mumps_icntl_8","ICNTL(8): scaling strategy (-2 to 8 or 77)","None",mumps->id.ICNTL(8),&mumps->id.ICNTL(8),NULL);
1242:   /* 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() */
1243:   PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",mumps->id.ICNTL(10),&mumps->id.ICNTL(10),NULL);
1244:   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);
1245:   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);
1246:   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);
1247:   PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage increase in the estimated working space","None",mumps->id.ICNTL(14),&mumps->id.ICNTL(14),NULL);
1248:   PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): computes the Schur complement","None",mumps->id.ICNTL(19),&mumps->id.ICNTL(19),NULL);
1249:   if (mumps->id.ICNTL(19) <= 0 || mumps->id.ICNTL(19) > 3) { /* reset any schur data (if any) */
1250:     MatMumpsResetSchur_Private(mumps);
1251:   }
1252:   /* 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 */
1253:   /* 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 */

1255:   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);
1256:   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);
1257:   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);
1258:   if (mumps->id.ICNTL(24)) {
1259:     mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */
1260:   }

1262:   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);
1263:   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);
1264:   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);
1265:   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);
1266:   PetscOptionsInt("-mat_mumps_icntl_29","ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis","None",mumps->id.ICNTL(29),&mumps->id.ICNTL(29),NULL);
1267:   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);
1268:   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);
1269:   /* 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 */
1270:   PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),NULL);

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

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

1280:   PetscOptionsIntArray("-mat_mumps_view_info","request INFO local to each processor","",info,&ninfo,NULL);
1281:   if (ninfo) {
1282:     if (ninfo > 40) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"number of INFO %d must <= 40\n",ninfo);
1283:     PetscMalloc1(ninfo,&mumps->info);
1284:     mumps->ninfo = ninfo;
1285:     for (i=0; i<ninfo; i++) {
1286:       if (info[i] < 0 || info[i]>40) {
1287:         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index of INFO %d must between 1 and 40\n",ninfo);
1288:       } else {
1289:         mumps->info[i] = info[i];
1290:       }
1291:     }
1292:   }

1294:   PetscOptionsEnd();
1295:   return(0);
1296: }

1300: PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps)
1301: {

1305:   MPI_Comm_rank(PetscObjectComm((PetscObject)A), &mumps->myid);
1306:   MPI_Comm_size(PetscObjectComm((PetscObject)A),&mumps->size);
1307:   MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(mumps->comm_mumps));

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

1311:   mumps->id.job = JOB_INIT;
1312:   mumps->id.par = 1;  /* host participates factorizaton and solve */
1313:   mumps->id.sym = mumps->sym;
1314:   PetscMUMPS_c(&mumps->id);

1316:   mumps->CleanUpMUMPS = PETSC_FALSE;
1317:   mumps->scat_rhs     = NULL;
1318:   mumps->scat_sol     = NULL;

1320:   /* set PETSc-MUMPS default options - override MUMPS default */
1321:   mumps->id.ICNTL(3) = 0;
1322:   mumps->id.ICNTL(4) = 0;
1323:   if (mumps->size == 1) {
1324:     mumps->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
1325:   } else {
1326:     mumps->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
1327:     mumps->id.ICNTL(20) = 0;   /* rhs is in dense format */
1328:     mumps->id.ICNTL(21) = 1;   /* distributed solution */
1329:   }

1331:   /* schur */
1332:   mumps->id.size_schur      = 0;
1333:   mumps->id.listvar_schur   = NULL;
1334:   mumps->id.schur           = NULL;
1335:   mumps->schur_second_solve = PETSC_FALSE;
1336:   mumps->sizeredrhs         = 0;
1337:   mumps->schur_pivots       = NULL;
1338:   mumps->schur_work         = NULL;
1339:   mumps->schur_sol          = NULL;
1340:   mumps->schur_sizesol      = 0;
1341:   mumps->schur_restored     = PETSC_TRUE;
1342:   mumps->schur_factored     = PETSC_FALSE;
1343:   mumps->schur_inverted     = PETSC_FALSE;
1344:   return(0);
1345: }

1347: /* Note Petsc r(=c) permutation is used when mumps->id.ICNTL(7)==1 with centralized assembled matrix input; otherwise r and c are ignored */
1350: PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1351: {
1352:   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
1354:   Vec            b;
1355:   IS             is_iden;
1356:   const PetscInt M = A->rmap->N;

1359:   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;

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

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

1366:   /* analysis phase */
1367:   /*----------------*/
1368:   mumps->id.job = JOB_FACTSYMBOLIC;
1369:   mumps->id.n   = M;
1370:   switch (mumps->id.ICNTL(18)) {
1371:   case 0:  /* centralized assembled matrix input */
1372:     if (!mumps->myid) {
1373:       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1374:       if (mumps->id.ICNTL(6)>1) {
1375:         mumps->id.a = (MumpsScalar*)mumps->val;
1376:       }
1377:       if (mumps->id.ICNTL(7) == 1) { /* use user-provide matrix ordering - assuming r = c ordering */
1378:         /*
1379:         PetscBool      flag;
1380:         ISEqual(r,c,&flag);
1381:         if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"row_perm != col_perm");
1382:         ISView(r,PETSC_VIEWER_STDOUT_SELF);
1383:          */
1384:         if (!mumps->myid) {
1385:           const PetscInt *idx;
1386:           PetscInt       i,*perm_in;

1388:           PetscMalloc1(M,&perm_in);
1389:           ISGetIndices(r,&idx);

1391:           mumps->id.perm_in = perm_in;
1392:           for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */
1393:           ISRestoreIndices(r,&idx);
1394:         }
1395:       }
1396:     }
1397:     break;
1398:   case 3:  /* distributed assembled matrix input (size>1) */
1399:     mumps->id.nz_loc = mumps->nz;
1400:     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1401:     if (mumps->id.ICNTL(6)>1) {
1402:       mumps->id.a_loc = (MumpsScalar*)mumps->val;
1403:     }
1404:     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1405:     if (!mumps->myid) {
1406:       VecCreateSeq(PETSC_COMM_SELF,A->rmap->N,&mumps->b_seq);
1407:       ISCreateStride(PETSC_COMM_SELF,A->rmap->N,0,1,&is_iden);
1408:     } else {
1409:       VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);
1410:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);
1411:     }
1412:     MatCreateVecs(A,NULL,&b);
1413:     VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);
1414:     ISDestroy(&is_iden);
1415:     VecDestroy(&b);
1416:     break;
1417:   }
1418:   PetscMUMPS_c(&mumps->id);
1419:   if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1));

1421:   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1422:   F->ops->solve           = MatSolve_MUMPS;
1423:   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1424:   F->ops->matsolve        = MatMatSolve_MUMPS;
1425:   return(0);
1426: }

1428: /* Note the Petsc r and c permutations are ignored */
1431: PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1432: {
1433:   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
1435:   Vec            b;
1436:   IS             is_iden;
1437:   const PetscInt M = A->rmap->N;

1440:   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;

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

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

1447:   /* analysis phase */
1448:   /*----------------*/
1449:   mumps->id.job = JOB_FACTSYMBOLIC;
1450:   mumps->id.n   = M;
1451:   switch (mumps->id.ICNTL(18)) {
1452:   case 0:  /* centralized assembled matrix input */
1453:     if (!mumps->myid) {
1454:       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1455:       if (mumps->id.ICNTL(6)>1) {
1456:         mumps->id.a = (MumpsScalar*)mumps->val;
1457:       }
1458:     }
1459:     break;
1460:   case 3:  /* distributed assembled matrix input (size>1) */
1461:     mumps->id.nz_loc = mumps->nz;
1462:     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1463:     if (mumps->id.ICNTL(6)>1) {
1464:       mumps->id.a_loc = (MumpsScalar*)mumps->val;
1465:     }
1466:     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1467:     if (!mumps->myid) {
1468:       VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);
1469:       ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);
1470:     } else {
1471:       VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);
1472:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);
1473:     }
1474:     MatCreateVecs(A,NULL,&b);
1475:     VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);
1476:     ISDestroy(&is_iden);
1477:     VecDestroy(&b);
1478:     break;
1479:   }
1480:   PetscMUMPS_c(&mumps->id);
1481:   if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1));

1483:   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1484:   F->ops->solve           = MatSolve_MUMPS;
1485:   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1486:   return(0);
1487: }

1489: /* Note the Petsc r permutation and factor info are ignored */
1492: PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
1493: {
1494:   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
1496:   Vec            b;
1497:   IS             is_iden;
1498:   const PetscInt M = A->rmap->N;

1501:   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;

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

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

1508:   /* analysis phase */
1509:   /*----------------*/
1510:   mumps->id.job = JOB_FACTSYMBOLIC;
1511:   mumps->id.n   = M;
1512:   switch (mumps->id.ICNTL(18)) {
1513:   case 0:  /* centralized assembled matrix input */
1514:     if (!mumps->myid) {
1515:       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1516:       if (mumps->id.ICNTL(6)>1) {
1517:         mumps->id.a = (MumpsScalar*)mumps->val;
1518:       }
1519:     }
1520:     break;
1521:   case 3:  /* distributed assembled matrix input (size>1) */
1522:     mumps->id.nz_loc = mumps->nz;
1523:     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1524:     if (mumps->id.ICNTL(6)>1) {
1525:       mumps->id.a_loc = (MumpsScalar*)mumps->val;
1526:     }
1527:     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1528:     if (!mumps->myid) {
1529:       VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);
1530:       ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);
1531:     } else {
1532:       VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);
1533:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);
1534:     }
1535:     MatCreateVecs(A,NULL,&b);
1536:     VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);
1537:     ISDestroy(&is_iden);
1538:     VecDestroy(&b);
1539:     break;
1540:   }
1541:   PetscMUMPS_c(&mumps->id);
1542:   if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1));

1544:   F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
1545:   F->ops->solve                 = MatSolve_MUMPS;
1546:   F->ops->solvetranspose        = MatSolve_MUMPS;
1547:   F->ops->matsolve              = MatMatSolve_MUMPS;
1548: #if defined(PETSC_USE_COMPLEX)
1549:   F->ops->getinertia = NULL;
1550: #else
1551:   F->ops->getinertia = MatGetInertia_SBAIJMUMPS;
1552: #endif
1553:   return(0);
1554: }

1558: PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
1559: {
1560:   PetscErrorCode    ierr;
1561:   PetscBool         iascii;
1562:   PetscViewerFormat format;
1563:   Mat_MUMPS         *mumps=(Mat_MUMPS*)A->spptr;

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

1569:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1570:   if (iascii) {
1571:     PetscViewerGetFormat(viewer,&format);
1572:     if (format == PETSC_VIEWER_ASCII_INFO) {
1573:       PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");
1574:       PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                   %d \n",mumps->id.sym);
1575:       PetscViewerASCIIPrintf(viewer,"  PAR (host participation):            %d \n",mumps->id.par);
1576:       PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):         %d \n",mumps->id.ICNTL(1));
1577:       PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg): %d \n",mumps->id.ICNTL(2));
1578:       PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):   %d \n",mumps->id.ICNTL(3));
1579:       PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):        %d \n",mumps->id.ICNTL(4));
1580:       PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):         %d \n",mumps->id.ICNTL(5));
1581:       PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):        %d \n",mumps->id.ICNTL(6));
1582:       PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (sequentia matrix ordering):%d \n",mumps->id.ICNTL(7));
1583:       PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scalling strategy):        %d \n",mumps->id.ICNTL(8));
1584:       PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements):  %d \n",mumps->id.ICNTL(10));
1585:       PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):          %d \n",mumps->id.ICNTL(11));
1586:       if (mumps->id.ICNTL(11)>0) {
1587:         PetscViewerASCIIPrintf(viewer,"    RINFOG(4) (inf norm of input mat):        %g\n",mumps->id.RINFOG(4));
1588:         PetscViewerASCIIPrintf(viewer,"    RINFOG(5) (inf norm of solution):         %g\n",mumps->id.RINFOG(5));
1589:         PetscViewerASCIIPrintf(viewer,"    RINFOG(6) (inf norm of residual):         %g\n",mumps->id.RINFOG(6));
1590:         PetscViewerASCIIPrintf(viewer,"    RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8));
1591:         PetscViewerASCIIPrintf(viewer,"    RINFOG(9) (error estimate):               %g \n",mumps->id.RINFOG(9));
1592:         PetscViewerASCIIPrintf(viewer,"    RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11));
1593:       }
1594:       PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",mumps->id.ICNTL(12));
1595:       PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",mumps->id.ICNTL(13));
1596:       PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",mumps->id.ICNTL(14));
1597:       /* ICNTL(15-17) not used */
1598:       PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",mumps->id.ICNTL(18));
1599:       PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Shur complement info):                       %d \n",mumps->id.ICNTL(19));
1600:       PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",mumps->id.ICNTL(20));
1601:       PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d \n",mumps->id.ICNTL(21));
1602:       PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",mumps->id.ICNTL(22));
1603:       PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",mumps->id.ICNTL(23));

1605:       PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",mumps->id.ICNTL(24));
1606:       PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",mumps->id.ICNTL(25));
1607:       PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",mumps->id.ICNTL(26));
1608:       PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",mumps->id.ICNTL(27));
1609:       PetscViewerASCIIPrintf(viewer,"  ICNTL(28) (use parallel or sequential ordering):        %d \n",mumps->id.ICNTL(28));
1610:       PetscViewerASCIIPrintf(viewer,"  ICNTL(29) (parallel ordering):                          %d \n",mumps->id.ICNTL(29));

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

1616:       PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",mumps->id.CNTL(1));
1617:       PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2));
1618:       PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",mumps->id.CNTL(3));
1619:       PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",mumps->id.CNTL(4));
1620:       PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",mumps->id.CNTL(5));

1622:       /* infomation local to each processor */
1623:       PetscViewerASCIIPrintf(viewer, "  RINFO(1) (local estimated flops for the elimination after analysis): \n");
1624:       PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);
1625:       PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %g \n",mumps->myid,mumps->id.RINFO(1));
1626:       PetscViewerFlush(viewer);
1627:       PetscViewerASCIIPrintf(viewer, "  RINFO(2) (local estimated flops for the assembly after factorization): \n");
1628:       PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(2));
1629:       PetscViewerFlush(viewer);
1630:       PetscViewerASCIIPrintf(viewer, "  RINFO(3) (local estimated flops for the elimination after factorization): \n");
1631:       PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(3));
1632:       PetscViewerFlush(viewer);

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

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

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

1646:       if (mumps->ninfo && mumps->ninfo <= 40){
1647:         PetscInt i;
1648:         for (i=0; i<mumps->ninfo; i++){
1649:           PetscViewerASCIIPrintf(viewer, "  INFO(%d): \n",mumps->info[i]);
1650:           PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(mumps->info[i]));
1651:           PetscViewerFlush(viewer);
1652:         }
1653:       }


1656:       PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);

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

1664:         PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(3));
1665:         PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(4));
1666:         PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",mumps->id.INFOG(5));
1667:         PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",mumps->id.INFOG(6));
1668:         PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively use after analysis): %d \n",mumps->id.INFOG(7));
1669:         PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",mumps->id.INFOG(8));
1670:         PetscViewerASCIIPrintf(viewer,"  INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",mumps->id.INFOG(9));
1671:         PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",mumps->id.INFOG(10));
1672:         PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",mumps->id.INFOG(11));
1673:         PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",mumps->id.INFOG(12));
1674:         PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",mumps->id.INFOG(13));
1675:         PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",mumps->id.INFOG(14));
1676:         PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",mumps->id.INFOG(15));
1677:         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));
1678:         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));
1679:         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));
1680:         PetscViewerASCIIPrintf(viewer,"  INFOG(19) (size of all MUMPS internal data allocated during factorization: sum over all processors): %d \n",mumps->id.INFOG(19));
1681:         PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",mumps->id.INFOG(20));
1682:         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));
1683:         PetscViewerASCIIPrintf(viewer,"  INFOG(22) (size in MB of memory effectively used during factorization - sum over all processors): %d \n",mumps->id.INFOG(22));
1684:         PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",mumps->id.INFOG(23));
1685:         PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",mumps->id.INFOG(24));
1686:         PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",mumps->id.INFOG(25));
1687:         PetscViewerASCIIPrintf(viewer,"  INFOG(28) (after factorization: number of null pivots encountered): %d\n",mumps->id.INFOG(28));
1688:         PetscViewerASCIIPrintf(viewer,"  INFOG(29) (after factorization: effective number of entries in the factors (sum over all processors)): %d\n",mumps->id.INFOG(29));
1689:         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));
1690:         PetscViewerASCIIPrintf(viewer,"  INFOG(32) (after analysis: type of analysis done): %d\n",mumps->id.INFOG(32));
1691:         PetscViewerASCIIPrintf(viewer,"  INFOG(33) (value used for ICNTL(8)): %d\n",mumps->id.INFOG(33));
1692:         PetscViewerASCIIPrintf(viewer,"  INFOG(34) (exponent of the determinant if determinant is requested): %d\n",mumps->id.INFOG(34));
1693:       }
1694:     }
1695:   }
1696:   return(0);
1697: }

1701: PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
1702: {
1703:   Mat_MUMPS *mumps =(Mat_MUMPS*)A->spptr;

1706:   info->block_size        = 1.0;
1707:   info->nz_allocated      = mumps->id.INFOG(20);
1708:   info->nz_used           = mumps->id.INFOG(20);
1709:   info->nz_unneeded       = 0.0;
1710:   info->assemblies        = 0.0;
1711:   info->mallocs           = 0.0;
1712:   info->memory            = 0.0;
1713:   info->fill_ratio_given  = 0;
1714:   info->fill_ratio_needed = 0;
1715:   info->factor_mallocs    = 0;
1716:   return(0);
1717: }

1719: /* -------------------------------------------------------------------------------------------*/
1722: PetscErrorCode MatMumpsSetSchurIndices_MUMPS(Mat F,PetscInt size,PetscInt idxs[])
1723: {
1724:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;

1728:   if (mumps->size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MUMPS parallel Schur complements not yet supported from PETSc\n");
1729:   if (mumps->id.size_schur != size) {
1730:     PetscFree2(mumps->id.listvar_schur,mumps->id.schur);
1731:     mumps->id.size_schur = size;
1732:     mumps->id.schur_lld = size;
1733:     PetscMalloc2(size,&mumps->id.listvar_schur,size*size,&mumps->id.schur);
1734:   }
1735:   PetscMemcpy(mumps->id.listvar_schur,idxs,size*sizeof(PetscInt));
1736:   if (F->factortype == MAT_FACTOR_LU) {
1737:     mumps->id.ICNTL(19) = 3; /* MUMPS returns full matrix */
1738:   } else {
1739:     mumps->id.ICNTL(19) = 2; /* MUMPS returns lower triangular part */
1740:   }
1741:   /* set a special value of ICNTL (not handled my MUMPS) to be used in the solve phase by PETSc */
1742:   mumps->id.ICNTL(26) = -1;
1743:   return(0);
1744: }

1748: /*@
1749:   MatMumpsSetSchurIndices - Set indices defining the Schur complement that MUMPS will compute during the factorization steps

1751:    Logically Collective on Mat

1753:    Input Parameters:
1754: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1755: .  size - size of the Schur complement indices
1756: -  idxs[] - array of Schur complement indices

1758:    Notes:
1759:    The user has to free the array idxs[] since the indices are copied by the routine.
1760:    MUMPS Schur complement mode is currently implemented for sequential matrices.

1762:    Level: advanced

1764:    References: MUMPS Users' Guide

1766: .seealso: MatGetFactor(), MatMumpsCreateSchurComplement(), MatMumpsGetSchurComplement()
1767: @*/
1768: PetscErrorCode MatMumpsSetSchurIndices(Mat F,PetscInt size,PetscInt idxs[])
1769: {

1775:   PetscTryMethod(F,"MatMumpsSetSchurIndices_C",(Mat,PetscInt,PetscInt[]),(F,size,idxs));
1776:   return(0);
1777: }

1779: /* -------------------------------------------------------------------------------------------*/
1782: PetscErrorCode MatMumpsCreateSchurComplement_MUMPS(Mat F,Mat* S)
1783: {
1784:   Mat            St;
1785:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
1786:   PetscScalar    *array;
1787: #if defined(PETSC_USE_COMPLEX)
1788:   PetscScalar    im = PetscSqrtScalar((PetscScalar)-1.0);
1789: #endif

1793:   if (!mumps->CleanUpMUMPS) { /* CleanUpMUMPS is set to true after numerical factorization */
1794:     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Numerical factorization phase not yet performed! You should call MatFactorSymbolic/Numeric before");
1795:   } else if (!mumps->id.ICNTL(19)) {
1796:     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatMumpsSetSchurIndices to enable it");
1797:   } else if (!mumps->id.size_schur) {
1798:     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatMumpsSetSchurIndices before");
1799:   } else if (!mumps->schur_restored) {
1800:     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur matrix has not been restored using MatMumpsRestoreSchurComplement");
1801:   }
1802:   MatCreate(PetscObjectComm((PetscObject)F),&St);
1803:   MatSetSizes(St,PETSC_DECIDE,PETSC_DECIDE,mumps->id.size_schur,mumps->id.size_schur);
1804:   MatSetType(St,MATDENSE);
1805:   MatSetUp(St);
1806:   MatDenseGetArray(St,&array);
1807:   if (!mumps->sym) { /* MUMPS always return a full matrix */
1808:     if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
1809:       PetscInt i,j,N=mumps->id.size_schur;
1810:       for (i=0;i<N;i++) {
1811:         for (j=0;j<N;j++) {
1812: #if !defined(PETSC_USE_COMPLEX)
1813:           PetscScalar val = mumps->id.schur[i*N+j];
1814: #else
1815:           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
1816: #endif
1817:           array[j*N+i] = val;
1818:         }
1819:       }
1820:     } else { /* stored by columns */
1821:       PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));
1822:     }
1823:   } else { /* either full or lower-triangular (not packed) */
1824:     if (mumps->id.ICNTL(19) == 2) { /* lower triangular stored by columns */
1825:       PetscInt i,j,N=mumps->id.size_schur;
1826:       for (i=0;i<N;i++) {
1827:         for (j=i;j<N;j++) {
1828: #if !defined(PETSC_USE_COMPLEX)
1829:           PetscScalar val = mumps->id.schur[i*N+j];
1830: #else
1831:           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
1832: #endif
1833:           array[i*N+j] = val;
1834:           array[j*N+i] = val;
1835:         }
1836:       }
1837:     } else if (mumps->id.ICNTL(19) == 3) { /* full matrix */
1838:       PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));
1839:     } else { /* ICNTL(19) == 1 lower triangular stored by rows */
1840:       PetscInt i,j,N=mumps->id.size_schur;
1841:       for (i=0;i<N;i++) {
1842:         for (j=0;j<i+1;j++) {
1843: #if !defined(PETSC_USE_COMPLEX)
1844:           PetscScalar val = mumps->id.schur[i*N+j];
1845: #else
1846:           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
1847: #endif
1848:           array[i*N+j] = val;
1849:           array[j*N+i] = val;
1850:         }
1851:       }
1852:     }
1853:   }
1854:   MatDenseRestoreArray(St,&array);
1855:   *S = St;
1856:   return(0);
1857: }

1861: /*@
1862:   MatMumpsCreateSchurComplement - Create a Schur complement matrix object using Schur data computed by MUMPS during the factorization step

1864:    Logically Collective on Mat

1866:    Input Parameters:
1867: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1868: .  *S - location where to return the Schur complement (MATDENSE)

1870:    Notes:
1871:    MUMPS Schur complement mode is currently implemented for sequential matrices.
1872:    The routine provides a copy of the Schur data stored within MUMPS data strutures. The caller must destroy the object when it is no longer needed.
1873:    If MatMumpsInvertSchurComplement has been called, the routine gets back the inverse

1875:    Level: advanced

1877:    References: MUMPS Users' Guide

1879: .seealso: MatGetFactor(), MatMumpsSetSchurIndices(), MatMumpsGetSchurComplement()
1880: @*/
1881: PetscErrorCode MatMumpsCreateSchurComplement(Mat F,Mat* S)
1882: {

1887:   PetscTryMethod(F,"MatMumpsCreateSchurComplement_C",(Mat,Mat*),(F,S));
1888:   return(0);
1889: }

1891: /* -------------------------------------------------------------------------------------------*/
1894: PetscErrorCode MatMumpsGetSchurComplement_MUMPS(Mat F,Mat* S)
1895: {
1896:   Mat            St;
1897:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;

1901:   if (!mumps->CleanUpMUMPS) { /* CleanUpMUMPS is set to true after numerical factorization */
1902:     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Numerical factorization phase not yet performed! You should call MatFactorSymbolic/Numeric before");
1903:   } else if (!mumps->id.ICNTL(19)) {
1904:     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatMumpsSetSchurIndices to enable it");
1905:   } else if (!mumps->id.size_schur) {
1906:     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatMumpsSetSchurIndices before");
1907:   } else if (!mumps->schur_restored) {
1908:     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur matrix has not been restored using MatMumpsRestoreSchurComplement");
1909:   }
1910:   /* It should be the responsibility of the user to handle different ICNTL(19) cases if they want to work with the raw data */
1911:   /* should I also add errors when the Schur complement has been already factored? */
1912:   MatCreateSeqDense(PetscObjectComm((PetscObject)F),mumps->id.size_schur,mumps->id.size_schur,(PetscScalar*)mumps->id.schur,&St);
1913:   *S = St;
1914:   mumps->schur_restored = PETSC_FALSE;
1915:   return(0);
1916: }

1920: /*@
1921:   MatMumpsGetSchurComplement - Get a Schur complement matrix object using the current status of the raw Schur data computed by MUMPS during the factorization step

1923:    Logically Collective on Mat

1925:    Input Parameters:
1926: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1927: .  *S - location where to return the Schur complement (MATDENSE)

1929:    Notes:
1930:    MUMPS Schur complement mode is currently implemented for sequential matrices.
1931:    The routine uses the pointer to the raw data of the Schur Complement stored within MUMPS data strutures. The caller should call MatMumpsRestoreSchurComplement when the object is no longer needed.
1932:    If MatMumpsInvertSchurComplement has been called, the routine gets back the inverse

1934:    Level: advanced

1936:    References: MUMPS Users' Guide

1938: .seealso: MatGetFactor(), MatMumpsSetSchurIndices(), MatMumpsRestoreSchurComplement(), MatMumpsCreateSchurComplement()
1939: @*/
1940: PetscErrorCode MatMumpsGetSchurComplement(Mat F,Mat* S)
1941: {

1946:   PetscUseMethod(F,"MatMumpsGetSchurComplement_C",(Mat,Mat*),(F,S));
1947:   return(0);
1948: }

1950: /* -------------------------------------------------------------------------------------------*/
1953: PetscErrorCode MatMumpsRestoreSchurComplement_MUMPS(Mat F,Mat* S)
1954: {
1955:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;

1959:   if (!mumps->CleanUpMUMPS) { /* CleanUpMUMPS is set to true after numerical factorization */
1960:     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Numerical factorization phase not yet performed! You should call MatFactorSymbolic/Numeric before");
1961:   } else if (!mumps->id.ICNTL(19)) {
1962:     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatMumpsSetSchurIndices to enable it");
1963:   } else if (!mumps->id.size_schur) {
1964:     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatMumpsSetSchurIndices before");
1965:   } else if (mumps->schur_restored) {
1966:     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur matrix has been already restored");
1967:   }
1968:   MatDestroy(S);
1969:   *S = NULL;
1970:   mumps->schur_restored = PETSC_TRUE;
1971:   return(0);
1972: }

1976: /*@
1977:   MatMumpsRestoreSchurComplement - Restore the Schur complement matrix object obtained from a call to MatGetSchurComplement

1979:    Logically Collective on Mat

1981:    Input Parameters:
1982: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1983: .  *S - location where the Schur complement is stored

1985:    Notes:
1986:    MUMPS Schur complement mode is currently implemented for sequential matrices.

1988:    Level: advanced

1990:    References: MUMPS Users' Guide

1992: .seealso: MatGetFactor(), MatMumpsSetSchurIndices(), MatMumpsGetSchurComplement(), MatMumpsCreateSchurComplement()
1993: @*/
1994: PetscErrorCode MatMumpsRestoreSchurComplement(Mat F,Mat* S)
1995: {

2000:   PetscUseMethod(F,"MatMumpsRestoreSchurComplement_C",(Mat,Mat*),(F,S));
2001:   return(0);
2002: }

2004: /* -------------------------------------------------------------------------------------------*/
2007: PetscErrorCode MatMumpsInvertSchurComplement_MUMPS(Mat F)
2008: {
2009:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;

2013:   if (!mumps->id.ICNTL(19)) { /* do nothing */
2014:     return(0);
2015:   }
2016:   if (!mumps->CleanUpMUMPS) { /* CleanUpMUMPS is set to true after numerical factorization */
2017:     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Numerical factorization phase not yet performed! You should call MatFactorSymbolic/Numeric before");
2018:   } else if (!mumps->id.size_schur) {
2019:     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatMumpsSetSchurIndices before");
2020:   } else if (!mumps->schur_restored) {
2021:     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur matrix has not been restored using MatMumpsRestoreSchurComplement");
2022:   }
2023:   MatMumpsInvertSchur_Private(mumps);
2024:   return(0);
2025: }

2029: /*@
2030:   MatMumpsInvertSchurComplement - Invert the raw Schur data computed by MUMPS during the factorization step

2032:    Logically Collective on Mat

2034:    Input Parameters:
2035: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface

2037:    Notes:
2038:    MUMPS Schur complement mode is currently implemented for sequential matrices.
2039:    The routine uses the pointer to the raw data of the Schur Complement stored within MUMPS data strutures.

2041:    Level: advanced

2043:    References: MUMPS Users' Guide

2045: .seealso: MatGetFactor(), MatMumpsSetSchurIndices()
2046: @*/
2047: PetscErrorCode MatMumpsInvertSchurComplement(Mat F)
2048: {

2053:   PetscTryMethod(F,"MatMumpsInvertSchurComplement_C",(Mat),(F));
2054:   return(0);
2055: }

2057: /* -------------------------------------------------------------------------------------------*/
2060: PetscErrorCode MatMumpsSolveSchurComplement_MUMPS(Mat F, Vec rhs, Vec sol)
2061: {
2062:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
2063:   MumpsScalar    *orhs;
2064:   PetscScalar    *osol,*nrhs,*nsol;
2065:   PetscInt       orhs_size,osol_size,olrhs_size;

2069:   if (!mumps->CleanUpMUMPS) { /* CleanUpMUMPS is set to true after numerical factorization */
2070:     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Numerical factorization phase not yet performed! You should call MatFactorSymbolic/Numeric before");
2071:   } else if (!mumps->id.ICNTL(19)) {
2072:     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatMumpsSetSchurIndices to enable it");
2073:   } else if (!mumps->id.size_schur) {
2074:     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatMumpsSetSchurIndices before");
2075:   } else if (!mumps->schur_restored) {
2076:     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur matrix has not been restored using MatMumpsRestoreSchurComplement");
2077:   }
2078:   /* swap pointers */
2079:   orhs = mumps->id.redrhs;
2080:   olrhs_size = mumps->id.lredrhs;
2081:   orhs_size = mumps->sizeredrhs;
2082:   osol = mumps->schur_sol;
2083:   osol_size = mumps->schur_sizesol;
2084:   VecGetArray(rhs,&nrhs);
2085:   VecGetArray(sol,&nsol);
2086:   mumps->id.redrhs = (MumpsScalar*)nrhs;
2087:   VecGetLocalSize(rhs,&mumps->sizeredrhs);
2088:   mumps->id.lredrhs = mumps->sizeredrhs;
2089:   mumps->schur_sol = nsol;
2090:   VecGetLocalSize(sol,&mumps->schur_sizesol);

2092:   /* solve Schur complement */
2093:   mumps->id.nrhs = 1;
2094:   MatMumpsSolveSchur_Private(mumps,PETSC_FALSE);
2095:   /* restore pointers */
2096:   VecRestoreArray(rhs,&nrhs);
2097:   VecRestoreArray(sol,&nsol);
2098:   mumps->id.redrhs = orhs;
2099:   mumps->id.lredrhs = olrhs_size;
2100:   mumps->sizeredrhs = orhs_size;
2101:   mumps->schur_sol = osol;
2102:   mumps->schur_sizesol = osol_size;
2103:   return(0);
2104: }

2108: /*@
2109:   MatMumpsSolveSchurComplement - Solve the Schur complement system computed by MUMPS during the factorization step

2111:    Logically Collective on Mat

2113:    Input Parameters:
2114: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2115: .  rhs - location where the right hand side of the Schur complement system is stored
2116: -  sol - location where the solution of the Schur complement system has to be returned

2118:    Notes:
2119:    MUMPS Schur complement mode is currently implemented for sequential matrices.
2120:    The sizes of the vectors should match the size of the Schur complement

2122:    Level: advanced

2124:    References: MUMPS Users' Guide

2126: .seealso: MatGetFactor(), MatMumpsSetSchurIndices()
2127: @*/
2128: PetscErrorCode MatMumpsSolveSchurComplement(Mat F, Vec rhs, Vec sol)
2129: {

2138:   PetscUseMethod(F,"MatMumpsSolveSchurComplement_C",(Mat,Vec,Vec),(F,rhs,sol));
2139:   return(0);
2140: }

2142: /* -------------------------------------------------------------------------------------------*/
2145: PetscErrorCode MatMumpsSolveSchurComplementTranspose_MUMPS(Mat F, Vec rhs, Vec sol)
2146: {
2147:   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
2148:   MumpsScalar    *orhs;
2149:   PetscScalar    *osol,*nrhs,*nsol;
2150:   PetscInt       orhs_size,osol_size;

2154:   if (!mumps->CleanUpMUMPS) { /* CleanUpMUMPS is set to true after numerical factorization */
2155:     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Numerical factorization phase not yet performed! You should call MatFactorSymbolic/Numeric before");
2156:   } else if (!mumps->id.ICNTL(19)) {
2157:     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatMumpsSetSchurIndices to enable it");
2158:   } else if (!mumps->id.size_schur) {
2159:     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatMumpsSetSchurIndices before");
2160:   } else if (!mumps->schur_restored) {
2161:     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur matrix has not been restored using MatMumpsRestoreSchurComplement");
2162:   }
2163:   /* swap pointers */
2164:   orhs = mumps->id.redrhs;
2165:   orhs_size = mumps->sizeredrhs;
2166:   osol = mumps->schur_sol;
2167:   osol_size = mumps->schur_sizesol;
2168:   VecGetArray(rhs,&nrhs);
2169:   VecGetArray(sol,&nsol);
2170:   mumps->id.redrhs = (MumpsScalar*)nrhs;
2171:   VecGetLocalSize(rhs,&mumps->sizeredrhs);
2172:   mumps->schur_sol = nsol;
2173:   VecGetLocalSize(sol,&mumps->schur_sizesol);

2175:   /* solve Schur complement */
2176:   mumps->id.nrhs = 1;
2177:   mumps->id.ICNTL(9) = 0;
2178:   MatMumpsSolveSchur_Private(mumps,PETSC_FALSE);
2179:   mumps->id.ICNTL(9) = 1;
2180:   /* restore pointers */
2181:   VecRestoreArray(rhs,&nrhs);
2182:   VecRestoreArray(sol,&nsol);
2183:   mumps->id.redrhs = orhs;
2184:   mumps->sizeredrhs = orhs_size;
2185:   mumps->schur_sol = osol;
2186:   mumps->schur_sizesol = osol_size;
2187:   return(0);
2188: }

2192: /*@
2193:   MatMumpsSolveSchurComplementTranspose - Solve the transpose of the Schur complement system computed by MUMPS during the factorization step

2195:    Logically Collective on Mat

2197:    Input Parameters:
2198: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2199: .  rhs - location where the right hand side of the Schur complement system is stored
2200: -  sol - location where the solution of the Schur complement system has to be returned

2202:    Notes:
2203:    MUMPS Schur complement mode is currently implemented for sequential matrices.
2204:    The sizes of the vectors should match the size of the Schur complement

2206:    Level: advanced

2208:    References: MUMPS Users' Guide

2210: .seealso: MatGetFactor(), MatMumpsSetSchurIndices()
2211: @*/
2212: PetscErrorCode MatMumpsSolveSchurComplementTranspose(Mat F, Vec rhs, Vec sol)
2213: {

2222:   PetscUseMethod(F,"MatMumpsSolveSchurComplementTranspose_C",(Mat,Vec,Vec),(F,rhs,sol));
2223:   return(0);
2224: }

2226: /* -------------------------------------------------------------------------------------------*/
2229: PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival)
2230: {
2231:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;

2234:   mumps->id.ICNTL(icntl) = ival;
2235:   return(0);
2236: }

2240: PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt *ival)
2241: {
2242:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;

2245:   *ival = mumps->id.ICNTL(icntl);
2246:   return(0);
2247: }

2251: /*@
2252:   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()

2254:    Logically Collective on Mat

2256:    Input Parameters:
2257: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2258: .  icntl - index of MUMPS parameter array ICNTL()
2259: -  ival - value of MUMPS ICNTL(icntl)

2261:   Options Database:
2262: .   -mat_mumps_icntl_<icntl> <ival>

2264:    Level: beginner

2266:    References: MUMPS Users' Guide

2268: .seealso: MatGetFactor()
2269: @*/
2270: PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival)
2271: {

2277:   PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));
2278:   return(0);
2279: }

2283: /*@
2284:   MatMumpsGetIcntl - Get MUMPS parameter ICNTL()

2286:    Logically Collective on Mat

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

2292:   Output Parameter:
2293: .  ival - value of MUMPS ICNTL(icntl)

2295:    Level: beginner

2297:    References: MUMPS Users' Guide

2299: .seealso: MatGetFactor()
2300: @*/
2301: PetscErrorCode MatMumpsGetIcntl(Mat F,PetscInt icntl,PetscInt *ival)
2302: {

2308:   PetscTryMethod(F,"MatMumpsGetIcntl_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));
2309:   return(0);
2310: }

2312: /* -------------------------------------------------------------------------------------------*/
2315: PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal val)
2316: {
2317:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;

2320:   mumps->id.CNTL(icntl) = val;
2321:   return(0);
2322: }

2326: PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal *val)
2327: {
2328:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;

2331:   *val = mumps->id.CNTL(icntl);
2332:   return(0);
2333: }

2337: /*@
2338:   MatMumpsSetCntl - Set MUMPS parameter CNTL()

2340:    Logically Collective on Mat

2342:    Input Parameters:
2343: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2344: .  icntl - index of MUMPS parameter array CNTL()
2345: -  val - value of MUMPS CNTL(icntl)

2347:   Options Database:
2348: .   -mat_mumps_cntl_<icntl> <val>

2350:    Level: beginner

2352:    References: MUMPS Users' Guide

2354: .seealso: MatGetFactor()
2355: @*/
2356: PetscErrorCode MatMumpsSetCntl(Mat F,PetscInt icntl,PetscReal val)
2357: {

2363:   PetscTryMethod(F,"MatMumpsSetCntl_C",(Mat,PetscInt,PetscReal),(F,icntl,val));
2364:   return(0);
2365: }

2369: /*@
2370:   MatMumpsGetCntl - Get MUMPS parameter CNTL()

2372:    Logically Collective on Mat

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

2378:   Output Parameter:
2379: .  val - value of MUMPS CNTL(icntl)

2381:    Level: beginner

2383:    References: MUMPS Users' Guide

2385: .seealso: MatGetFactor()
2386: @*/
2387: PetscErrorCode MatMumpsGetCntl(Mat F,PetscInt icntl,PetscReal *val)
2388: {

2394:   PetscTryMethod(F,"MatMumpsGetCntl_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));
2395:   return(0);
2396: }

2400: PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F,PetscInt icntl,PetscInt *info)
2401: {
2402:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;

2405:   *info = mumps->id.INFO(icntl);
2406:   return(0);
2407: }

2411: PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F,PetscInt icntl,PetscInt *infog)
2412: {
2413:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;

2416:   *infog = mumps->id.INFOG(icntl);
2417:   return(0);
2418: }

2422: PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfo)
2423: {
2424:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;

2427:   *rinfo = mumps->id.RINFO(icntl);
2428:   return(0);
2429: }

2433: PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfog)
2434: {
2435:   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;

2438:   *rinfog = mumps->id.RINFOG(icntl);
2439:   return(0);
2440: }

2444: /*@
2445:   MatMumpsGetInfo - Get MUMPS parameter INFO()

2447:    Logically Collective on Mat

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

2453:   Output Parameter:
2454: .  ival - value of MUMPS INFO(icntl)

2456:    Level: beginner

2458:    References: MUMPS Users' Guide

2460: .seealso: MatGetFactor()
2461: @*/
2462: PetscErrorCode MatMumpsGetInfo(Mat F,PetscInt icntl,PetscInt *ival)
2463: {

2468:   PetscTryMethod(F,"MatMumpsGetInfo_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));
2469:   return(0);
2470: }

2474: /*@
2475:   MatMumpsGetInfog - Get MUMPS parameter INFOG()

2477:    Logically Collective on Mat

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

2483:   Output Parameter:
2484: .  ival - value of MUMPS INFOG(icntl)

2486:    Level: beginner

2488:    References: MUMPS Users' Guide

2490: .seealso: MatGetFactor()
2491: @*/
2492: PetscErrorCode MatMumpsGetInfog(Mat F,PetscInt icntl,PetscInt *ival)
2493: {

2498:   PetscTryMethod(F,"MatMumpsGetInfog_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));
2499:   return(0);
2500: }

2504: /*@
2505:   MatMumpsGetRinfo - Get MUMPS parameter RINFO()

2507:    Logically Collective on Mat

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

2513:   Output Parameter:
2514: .  val - value of MUMPS RINFO(icntl)

2516:    Level: beginner

2518:    References: MUMPS Users' Guide

2520: .seealso: MatGetFactor()
2521: @*/
2522: PetscErrorCode MatMumpsGetRinfo(Mat F,PetscInt icntl,PetscReal *val)
2523: {

2528:   PetscTryMethod(F,"MatMumpsGetRinfo_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));
2529:   return(0);
2530: }

2534: /*@
2535:   MatMumpsGetRinfog - Get MUMPS parameter RINFOG()

2537:    Logically Collective on Mat

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

2543:   Output Parameter:
2544: .  val - value of MUMPS RINFOG(icntl)

2546:    Level: beginner

2548:    References: MUMPS Users' Guide

2550: .seealso: MatGetFactor()
2551: @*/
2552: PetscErrorCode MatMumpsGetRinfog(Mat F,PetscInt icntl,PetscReal *val)
2553: {

2558:   PetscTryMethod(F,"MatMumpsGetRinfog_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));
2559:   return(0);
2560: }

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

2566:   Works with MATAIJ and MATSBAIJ matrices

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

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

2572:   Options Database Keys:
2573: +  -mat_mumps_icntl_1 <6>: ICNTL(1): output stream for error messages (None)
2574: .  -mat_mumps_icntl_2 <0>: ICNTL(2): output stream for diagnostic printing, statistics, and warning (None)
2575: .  -mat_mumps_icntl_3 <0>: ICNTL(3): output stream for global information, collected on the host (None)
2576: .  -mat_mumps_icntl_4 <0>: ICNTL(4): level of printing (0 to 4) (None)
2577: .  -mat_mumps_icntl_6 <7>: ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7) (None)
2578: .  -mat_mumps_icntl_7 <7>: ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis (None)
2579: .  -mat_mumps_icntl_8 <77>: ICNTL(8): scaling strategy (-2 to 8 or 77) (None)
2580: .  -mat_mumps_icntl_10 <0>: ICNTL(10): max num of refinements (None)
2581: .  -mat_mumps_icntl_11 <0>: ICNTL(11): statistics related to an error analysis (via -ksp_view) (None)
2582: .  -mat_mumps_icntl_12 <1>: ICNTL(12): an ordering strategy for symmetric matrices (0 to 3) (None)
2583: .  -mat_mumps_icntl_13 <0>: ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting (None)
2584: .  -mat_mumps_icntl_14 <20>: ICNTL(14): percentage increase in the estimated working space (None)
2585: .  -mat_mumps_icntl_19 <0>: ICNTL(19): computes the Schur complement (None)
2586: .  -mat_mumps_icntl_22 <0>: ICNTL(22): in-core/out-of-core factorization and solve (0 or 1) (None)
2587: .  -mat_mumps_icntl_23 <0>: ICNTL(23): max size of the working memory (MB) that can allocate per processor (None)
2588: .  -mat_mumps_icntl_24 <0>: ICNTL(24): detection of null pivot rows (0 or 1) (None)
2589: .  -mat_mumps_icntl_25 <0>: ICNTL(25): compute a solution of a deficient matrix and a null space basis (None)
2590: .  -mat_mumps_icntl_26 <0>: ICNTL(26): drives the solution phase if a Schur complement matrix (None)
2591: .  -mat_mumps_icntl_28 <1>: ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering (None)
2592: .  -mat_mumps_icntl_29 <0>: ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis (None)
2593: .  -mat_mumps_icntl_30 <0>: ICNTL(30): compute user-specified set of entries in inv(A) (None)
2594: .  -mat_mumps_icntl_31 <0>: ICNTL(31): indicates which factors may be discarded during factorization (None)
2595: .  -mat_mumps_icntl_33 <0>: ICNTL(33): compute determinant (None)
2596: .  -mat_mumps_cntl_1 <0.01>: CNTL(1): relative pivoting threshold (None)
2597: .  -mat_mumps_cntl_2 <1.49012e-08>: CNTL(2): stopping criterion of refinement (None)
2598: .  -mat_mumps_cntl_3 <0>: CNTL(3): absolute pivoting threshold (None)
2599: .  -mat_mumps_cntl_4 <-1>: CNTL(4): value for static pivoting (None)
2600: -  -mat_mumps_cntl_5 <0>: CNTL(5): fixation for null pivots (None)

2602:   Level: beginner

2604: .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage

2606: M*/

2610: static PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type)
2611: {
2613:   *type = MATSOLVERMUMPS;
2614:   return(0);
2615: }

2617: /* MatGetFactor for Seq and MPI AIJ matrices */
2620: PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
2621: {
2622:   Mat            B;
2624:   Mat_MUMPS      *mumps;
2625:   PetscBool      isSeqAIJ;

2628:   /* Create the factorization matrix */
2629:   PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);
2630:   MatCreate(PetscObjectComm((PetscObject)A),&B);
2631:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
2632:   MatSetType(B,((PetscObject)A)->type_name);
2633:   if (isSeqAIJ) {
2634:     MatSeqAIJSetPreallocation(B,0,NULL);
2635:   } else {
2636:     MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);
2637:   }

2639:   PetscNewLog(B,&mumps);

2641:   B->ops->view        = MatView_MUMPS;
2642:   B->ops->getinfo     = MatGetInfo_MUMPS;
2643:   B->ops->getdiagonal = MatGetDiagonal_MUMPS;

2645:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);
2646:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);
2647:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);
2648:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);
2649:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);

2651:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);
2652:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);
2653:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);
2654:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);

2656:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetSchurIndices_C",MatMumpsSetSchurIndices_MUMPS);
2657:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsInvertSchurComplement_C",MatMumpsInvertSchurComplement_MUMPS);
2658:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsCreateSchurComplement_C",MatMumpsCreateSchurComplement_MUMPS);
2659:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetSchurComplement_C",MatMumpsGetSchurComplement_MUMPS);
2660:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsRestoreSchurComplement_C",MatMumpsRestoreSchurComplement_MUMPS);
2661:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSolveSchurComplement_C",MatMumpsSolveSchurComplement_MUMPS);
2662:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSolveSchurComplementTranspose_C",MatMumpsSolveSchurComplementTranspose_MUMPS);

2664:   if (ftype == MAT_FACTOR_LU) {
2665:     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
2666:     B->factortype            = MAT_FACTOR_LU;
2667:     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
2668:     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
2669:     mumps->sym = 0;
2670:   } else {
2671:     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2672:     B->factortype                  = MAT_FACTOR_CHOLESKY;
2673:     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
2674:     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
2675: #if defined(PETSC_USE_COMPLEX)
2676:     mumps->sym = 2;
2677: #else
2678:     if (A->spd_set && A->spd) mumps->sym = 1;
2679:     else                      mumps->sym = 2;
2680: #endif
2681:   }

2683:   mumps->isAIJ    = PETSC_TRUE;
2684:   mumps->Destroy  = B->ops->destroy;
2685:   B->ops->destroy = MatDestroy_MUMPS;
2686:   B->spptr        = (void*)mumps;

2688:   PetscInitializeMUMPS(A,mumps);

2690:   *F = B;
2691:   return(0);
2692: }

2694: /* MatGetFactor for Seq and MPI SBAIJ matrices */
2697: PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
2698: {
2699:   Mat            B;
2701:   Mat_MUMPS      *mumps;
2702:   PetscBool      isSeqSBAIJ;

2705:   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
2706:   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");
2707:   PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);
2708:   /* Create the factorization matrix */
2709:   MatCreate(PetscObjectComm((PetscObject)A),&B);
2710:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
2711:   MatSetType(B,((PetscObject)A)->type_name);
2712:   PetscNewLog(B,&mumps);
2713:   if (isSeqSBAIJ) {
2714:     MatSeqSBAIJSetPreallocation(B,1,0,NULL);

2716:     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
2717:   } else {
2718:     MatMPISBAIJSetPreallocation(B,1,0,NULL,0,NULL);

2720:     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
2721:   }

2723:   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2724:   B->ops->view                   = MatView_MUMPS;
2725:   B->ops->getdiagonal            = MatGetDiagonal_MUMPS;

2727:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);
2728:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);
2729:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);
2730:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);
2731:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);

2733:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);
2734:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);
2735:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);
2736:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);

2738:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetSchurIndices_C",MatMumpsSetSchurIndices_MUMPS);
2739:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsInvertSchurComplement_C",MatMumpsInvertSchurComplement_MUMPS);
2740:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsCreateSchurComplement_C",MatMumpsCreateSchurComplement_MUMPS);
2741:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetSchurComplement_C",MatMumpsGetSchurComplement_MUMPS);
2742:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsRestoreSchurComplement_C",MatMumpsRestoreSchurComplement_MUMPS);
2743:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSolveSchurComplement_C",MatMumpsSolveSchurComplement_MUMPS);
2744:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSolveSchurComplementTranspose_C",MatMumpsSolveSchurComplementTranspose_MUMPS);

2746:   B->factortype = MAT_FACTOR_CHOLESKY;
2747: #if defined(PETSC_USE_COMPLEX)
2748:   mumps->sym = 2;
2749: #else
2750:   if (A->spd_set && A->spd) mumps->sym = 1;
2751:   else                      mumps->sym = 2;
2752: #endif

2754:   mumps->isAIJ    = PETSC_FALSE;
2755:   mumps->Destroy  = B->ops->destroy;
2756:   B->ops->destroy = MatDestroy_MUMPS;
2757:   B->spptr        = (void*)mumps;

2759:   PetscInitializeMUMPS(A,mumps);

2761:   *F = B;
2762:   return(0);
2763: }

2767: PETSC_EXTERN PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
2768: {
2769:   Mat            B;
2771:   Mat_MUMPS      *mumps;
2772:   PetscBool      isSeqBAIJ;

2775:   /* Create the factorization matrix */
2776:   PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);
2777:   MatCreate(PetscObjectComm((PetscObject)A),&B);
2778:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
2779:   MatSetType(B,((PetscObject)A)->type_name);
2780:   if (isSeqBAIJ) {
2781:     MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,NULL);
2782:   } else {
2783:     MatMPIBAIJSetPreallocation(B,A->rmap->bs,0,NULL,0,NULL);
2784:   }

2786:   PetscNewLog(B,&mumps);
2787:   if (ftype == MAT_FACTOR_LU) {
2788:     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
2789:     B->factortype            = MAT_FACTOR_LU;
2790:     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
2791:     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
2792:     mumps->sym = 0;
2793:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n");

2795:   B->ops->view        = MatView_MUMPS;
2796:   B->ops->getdiagonal = MatGetDiagonal_MUMPS;

2798:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);
2799:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);
2800:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);
2801:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);
2802:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);

2804:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);
2805:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);
2806:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);
2807:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);

2809:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetSchurIndices_C",MatMumpsSetSchurIndices_MUMPS);
2810:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsInvertSchurComplement_C",MatMumpsInvertSchurComplement_MUMPS);
2811:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsCreateSchurComplement_C",MatMumpsCreateSchurComplement_MUMPS);
2812:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetSchurComplement_C",MatMumpsGetSchurComplement_MUMPS);
2813:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsRestoreSchurComplement_C",MatMumpsRestoreSchurComplement_MUMPS);
2814:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSolveSchurComplement_C",MatMumpsSolveSchurComplement_MUMPS);
2815:   PetscObjectComposeFunction((PetscObject)B,"MatMumpsSolveSchurComplementTranspose_C",MatMumpsSolveSchurComplementTranspose_MUMPS);

2817:   mumps->isAIJ    = PETSC_TRUE;
2818:   mumps->Destroy  = B->ops->destroy;
2819:   B->ops->destroy = MatDestroy_MUMPS;
2820:   B->spptr        = (void*)mumps;

2822:   PetscInitializeMUMPS(A,mumps);

2824:   *F = B;
2825:   return(0);
2826: }

2828: PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mumps(Mat,MatFactorType,Mat*);
2829: PETSC_EXTERN PetscErrorCode MatGetFactor_baij_mumps(Mat,MatFactorType,Mat*);
2830: PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat,MatFactorType,Mat*);

2834: PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_MUMPS(void)
2835: {

2839:   MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);
2840:   MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);
2841:   MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);
2842:   MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);
2843:   MatSolverPackageRegister(MATSOLVERMUMPS,MATMPISBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);
2844:   MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);
2845:   MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);
2846:   MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);
2847:   MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);
2848:   MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);
2849:   return(0);
2850: }