Actual source code: mumps.c

petsc-3.3-p7 2013-05-11
  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>

  9: EXTERN_C_BEGIN
 10: #if defined(PETSC_USE_COMPLEX)
 11: #include <zmumps_c.h>
 12: #else
 13: #include <dmumps_c.h> 
 14: #endif
 15: EXTERN_C_END
 16: #define JOB_INIT -1
 17: #define JOB_FACTSYMBOLIC 1
 18: #define JOB_FACTNUMERIC 2
 19: #define JOB_SOLVE 3
 20: #define JOB_END -2


 23: /* macros s.t. indices match MUMPS documentation */
 24: #define ICNTL(I) icntl[(I)-1]
 25: #define CNTL(I) cntl[(I)-1] 
 26: #define INFOG(I) infog[(I)-1]
 27: #define INFO(I) info[(I)-1]
 28: #define RINFOG(I) rinfog[(I)-1]
 29: #define RINFO(I) rinfo[(I)-1]

 31: typedef struct {
 32: #if defined(PETSC_USE_COMPLEX)
 33:   ZMUMPS_STRUC_C id;
 34: #else
 35:   DMUMPS_STRUC_C id;
 36: #endif
 37:   MatStructure   matstruc;
 38:   PetscMPIInt    myid,size;
 39:   PetscInt       *irn,*jcn,nz,sym,nSolve;
 40:   PetscScalar    *val;
 41:   MPI_Comm       comm_mumps;
 42:   VecScatter     scat_rhs, scat_sol;
 43:   PetscBool      isAIJ,CleanUpMUMPS;
 44:   Vec            b_seq,x_seq;
 45:   PetscErrorCode (*Destroy)(Mat);
 46:   PetscErrorCode (*ConvertToTriples)(Mat, int, MatReuse, int*, int**, int**, PetscScalar**);
 47: } Mat_MUMPS;

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


 52: /* MatConvertToTriples_A_B */
 53: /*convert Petsc matrix to triples: row[nz], col[nz], val[nz] */
 54: /*
 55:   input: 
 56:     A       - matrix in aij,baij or sbaij (bs=1) format
 57:     shift   - 0: C style output triple; 1: Fortran style output triple.
 58:     reuse   - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple  
 59:               MAT_REUSE_MATRIX:   only the values in v array are updated
 60:   output:     
 61:     nnz     - dim of r, c, and v (number of local nonzero entries of A)
 62:     r, c, v - row and col index, matrix values (matrix triples) 
 63:  */

 67: PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
 68: {
 69:   const PetscInt   *ai,*aj,*ajj,M=A->rmap->n;
 70:   PetscInt         nz,rnz,i,j;
 71:   PetscErrorCode   ierr;
 72:   PetscInt         *row,*col;
 73:   Mat_SeqAIJ       *aa=(Mat_SeqAIJ*)A->data;

 76:   *v=aa->a;
 77:   if (reuse == MAT_INITIAL_MATRIX){
 78:     nz = aa->nz; ai = aa->i; aj = aa->j;
 79:     *nnz = nz;
 80:     PetscMalloc(2*nz*sizeof(PetscInt), &row);
 81:     col  = row + nz;

 83:     nz = 0;
 84:     for(i=0; i<M; i++) {
 85:       rnz = ai[i+1] - ai[i];
 86:       ajj = aj + ai[i];
 87:       for(j=0; j<rnz; j++) {
 88:         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
 89:       }
 90:     }
 91:     *r = row; *c = col;
 92:   }
 93:   return(0);
 94: }

 98: PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
 99: {
100:   Mat_SeqBAIJ        *aa=(Mat_SeqBAIJ*)A->data;
101:   const PetscInt     *ai,*aj,*ajj,bs=A->rmap->bs,bs2=aa->bs2,M=A->rmap->N/bs;
102:   PetscInt           nz,idx=0,rnz,i,j,k,m;
103:   PetscErrorCode     ierr;
104:   PetscInt           *row,*col;

107:   *v = aa->a;
108:   if (reuse == MAT_INITIAL_MATRIX){
109:     ai = aa->i; aj = aa->j;
110:     nz = bs2*aa->nz;
111:     *nnz = nz;
112:     PetscMalloc(2*nz*sizeof(PetscInt), &row);
113:     col  = row + nz;

115:     for(i=0; i<M; i++) {
116:       ajj = aj + ai[i];
117:       rnz = ai[i+1] - ai[i];
118:       for(k=0; k<rnz; k++) {
119:         for(j=0; j<bs; j++) {
120:           for(m=0; m<bs; m++) {
121:             row[idx]     = i*bs + m + shift;
122:             col[idx++]   = bs*(ajj[k]) + j + shift;
123:           }
124:         }
125:       }
126:     }
127:     *r = row; *c = col;
128:   }
129:   return(0);
130: }

134: PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
135: {
136:   const PetscInt   *ai, *aj,*ajj,M=A->rmap->n;
137:   PetscInt         nz,rnz,i,j;
138:   PetscErrorCode   ierr;
139:   PetscInt         *row,*col;
140:   Mat_SeqSBAIJ     *aa=(Mat_SeqSBAIJ*)A->data;

143:   *v = aa->a;
144:   if (reuse == MAT_INITIAL_MATRIX){
145:     nz = aa->nz;ai=aa->i; aj=aa->j;*v=aa->a;
146:     *nnz = nz;
147:     PetscMalloc(2*nz*sizeof(PetscInt), &row);
148:     col  = row + nz;

150:     nz = 0;
151:     for(i=0; i<M; i++) {
152:       rnz = ai[i+1] - ai[i];
153:       ajj = aj + ai[i];
154:       for(j=0; j<rnz; j++) {
155:         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
156:       }
157:     }
158:     *r = row; *c = col;
159:   }
160:   return(0);
161: }

165: PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
166: {
167:   const PetscInt     *ai,*aj,*ajj,*adiag,M=A->rmap->n;
168:   PetscInt           nz,rnz,i,j;
169:   const PetscScalar  *av,*v1;
170:   PetscScalar        *val;
171:   PetscErrorCode     ierr;
172:   PetscInt           *row,*col;
173:   Mat_SeqSBAIJ       *aa=(Mat_SeqSBAIJ*)A->data;

176:   ai=aa->i; aj=aa->j;av=aa->a;
177:   adiag=aa->diag;
178:   if (reuse == MAT_INITIAL_MATRIX){
179:     nz = M + (aa->nz-M)/2;
180:     *nnz = nz;
181:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
182:     col  = row + nz;
183:     val  = (PetscScalar*)(col + nz);

185:     nz = 0;
186:     for(i=0; i<M; i++) {
187:       rnz = ai[i+1] - adiag[i];
188:       ajj  = aj + adiag[i];
189:       v1   = av + adiag[i];
190:       for(j=0; j<rnz; j++) {
191:         row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
192:       }
193:     }
194:     *r = row; *c = col; *v = val;
195:   } else {
196:     nz = 0; val = *v;
197:     for(i=0; i <M; i++) {
198:       rnz = ai[i+1] - adiag[i];
199:       ajj = aj + adiag[i];
200:       v1  = av + adiag[i];
201:       for(j=0; j<rnz; j++) {
202:         val[nz++] = v1[j];
203:       }
204:     }
205:   }
206:   return(0);
207: }

211: PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
212: {
213:   const PetscInt     *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
214:   PetscErrorCode     ierr;
215:   PetscInt           rstart,nz,i,j,jj,irow,countA,countB;
216:   PetscInt           *row,*col;
217:   const PetscScalar  *av, *bv,*v1,*v2;
218:   PetscScalar        *val;
219:   Mat_MPISBAIJ       *mat =  (Mat_MPISBAIJ*)A->data;
220:   Mat_SeqSBAIJ       *aa=(Mat_SeqSBAIJ*)(mat->A)->data;
221:   Mat_SeqBAIJ        *bb=(Mat_SeqBAIJ*)(mat->B)->data;

224:   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
225:   garray = mat->garray;
226:   av=aa->a; bv=bb->a;

228:   if (reuse == MAT_INITIAL_MATRIX){
229:     nz = aa->nz + bb->nz;
230:     *nnz = nz;
231:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
232:     col  = row + nz;
233:     val  = (PetscScalar*)(col + nz);

235:     *r = row; *c = col; *v = val;
236:   } else {
237:     row = *r; col = *c; val = *v;
238:   }

240:   jj = 0; irow = rstart;
241:   for ( i=0; i<m; i++ ) {
242:     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
243:     countA = ai[i+1] - ai[i];
244:     countB = bi[i+1] - bi[i];
245:     bjj    = bj + bi[i];
246:     v1     = av + ai[i];
247:     v2     = bv + bi[i];

249:     /* A-part */
250:     for (j=0; j<countA; j++){
251:       if (reuse == MAT_INITIAL_MATRIX) {
252:         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
253:       }
254:       val[jj++] = v1[j];
255:     }

257:     /* B-part */
258:     for(j=0; j < countB; j++){
259:       if (reuse == MAT_INITIAL_MATRIX) {
260:         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
261:       }
262:       val[jj++] = v2[j];
263:     }
264:     irow++;
265:   }
266:   return(0);
267: }

271: PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
272: {
273:   const PetscInt     *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
274:   PetscErrorCode     ierr;
275:   PetscInt           rstart,nz,i,j,jj,irow,countA,countB;
276:   PetscInt           *row,*col;
277:   const PetscScalar  *av, *bv,*v1,*v2;
278:   PetscScalar        *val;
279:   Mat_MPIAIJ         *mat =  (Mat_MPIAIJ*)A->data;
280:   Mat_SeqAIJ         *aa=(Mat_SeqAIJ*)(mat->A)->data;
281:   Mat_SeqAIJ         *bb=(Mat_SeqAIJ*)(mat->B)->data;

284:   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
285:   garray = mat->garray;
286:   av=aa->a; bv=bb->a;

288:   if (reuse == MAT_INITIAL_MATRIX){
289:     nz = aa->nz + bb->nz;
290:     *nnz = nz;
291:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
292:     col  = row + nz;
293:     val  = (PetscScalar*)(col + nz);

295:     *r = row; *c = col; *v = val;
296:   } else {
297:     row = *r; col = *c; val = *v;
298:   }

300:   jj = 0; irow = rstart;
301:   for ( i=0; i<m; i++ ) {
302:     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
303:     countA = ai[i+1] - ai[i];
304:     countB = bi[i+1] - bi[i];
305:     bjj    = bj + bi[i];
306:     v1     = av + ai[i];
307:     v2     = bv + bi[i];

309:     /* A-part */
310:     for (j=0; j<countA; j++){
311:       if (reuse == MAT_INITIAL_MATRIX){
312:         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
313:       }
314:       val[jj++] = v1[j];
315:     }

317:     /* B-part */
318:     for(j=0; j < countB; j++){
319:       if (reuse == MAT_INITIAL_MATRIX){
320:         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
321:       }
322:       val[jj++] = v2[j];
323:     }
324:     irow++;
325:   }
326:   return(0);
327: }

331: PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
332: {
333:   Mat_MPIBAIJ        *mat =  (Mat_MPIBAIJ*)A->data;
334:   Mat_SeqBAIJ        *aa=(Mat_SeqBAIJ*)(mat->A)->data;
335:   Mat_SeqBAIJ        *bb=(Mat_SeqBAIJ*)(mat->B)->data;
336:   const PetscInt     *ai = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj;
337:   const PetscInt     *garray = mat->garray,mbs=mat->mbs,rstart=A->rmap->rstart;
338:   const PetscInt     bs = A->rmap->bs,bs2=mat->bs2;
339:   PetscErrorCode     ierr;
340:   PetscInt           nz,i,j,k,n,jj,irow,countA,countB,idx;
341:   PetscInt           *row,*col;
342:   const PetscScalar  *av=aa->a, *bv=bb->a,*v1,*v2;
343:   PetscScalar        *val;


347:   if (reuse == MAT_INITIAL_MATRIX) {
348:     nz = bs2*(aa->nz + bb->nz);
349:     *nnz = nz;
350:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
351:     col  = row + nz;
352:     val  = (PetscScalar*)(col + nz);

354:     *r = row; *c = col; *v = val;
355:   } else {
356:     row = *r; col = *c; val = *v;
357:   }

359:   jj = 0; irow = rstart;
360:   for ( i=0; i<mbs; i++ ) {
361:     countA = ai[i+1] - ai[i];
362:     countB = bi[i+1] - bi[i];
363:     ajj    = aj + ai[i];
364:     bjj    = bj + bi[i];
365:     v1     = av + bs2*ai[i];
366:     v2     = bv + bs2*bi[i];

368:     idx = 0;
369:     /* A-part */
370:     for (k=0; k<countA; k++){
371:       for (j=0; j<bs; j++) {
372:         for (n=0; n<bs; n++) {
373:           if (reuse == MAT_INITIAL_MATRIX){
374:             row[jj] = irow + n + shift;
375:             col[jj] = rstart + bs*ajj[k] + j + shift;
376:           }
377:           val[jj++] = v1[idx++];
378:         }
379:       }
380:     }

382:     idx = 0;
383:     /* B-part */
384:     for(k=0; k<countB; k++){
385:       for (j=0; j<bs; j++) {
386:         for (n=0; n<bs; n++) {
387:           if (reuse == MAT_INITIAL_MATRIX){
388:             row[jj] = irow + n + shift;
389:             col[jj] = bs*garray[bjj[k]] + j + shift;
390:           }
391:           val[jj++] = v2[idx++];
392:         }
393:       }
394:     }
395:     irow += bs;
396:   }
397:   return(0);
398: }

402: PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
403: {
404:   const PetscInt     *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
405:   PetscErrorCode     ierr;
406:   PetscInt           rstart,nz,nza,nzb,i,j,jj,irow,countA,countB;
407:   PetscInt           *row,*col;
408:   const PetscScalar  *av, *bv,*v1,*v2;
409:   PetscScalar        *val;
410:   Mat_MPIAIJ         *mat =  (Mat_MPIAIJ*)A->data;
411:   Mat_SeqAIJ         *aa=(Mat_SeqAIJ*)(mat->A)->data;
412:   Mat_SeqAIJ         *bb=(Mat_SeqAIJ*)(mat->B)->data;

415:   ai=aa->i; aj=aa->j; adiag=aa->diag;
416:   bi=bb->i; bj=bb->j; garray = mat->garray;
417:   av=aa->a; bv=bb->a;
418:   rstart = A->rmap->rstart;

420:   if (reuse == MAT_INITIAL_MATRIX) {
421:     nza = 0;    /* num of upper triangular entries in mat->A, including diagonals */
422:     nzb = 0;    /* num of upper triangular entries in mat->B */
423:     for(i=0; i<m; i++){
424:       nza    += (ai[i+1] - adiag[i]);
425:       countB  = bi[i+1] - bi[i];
426:       bjj     = bj + bi[i];
427:       for (j=0; j<countB; j++){
428:         if (garray[bjj[j]] > rstart) nzb++;
429:       }
430:     }
431: 
432:     nz = nza + nzb; /* total nz of upper triangular part of mat */
433:     *nnz = nz;
434:     PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
435:     col  = row + nz;
436:     val  = (PetscScalar*)(col + nz);

438:     *r = row; *c = col; *v = val;
439:   } else {
440:     row = *r; col = *c; val = *v;
441:   }

443:   jj = 0; irow = rstart;
444:   for ( i=0; i<m; i++ ) {
445:     ajj    = aj + adiag[i];                 /* ptr to the beginning of the diagonal of this row */
446:     v1     = av + adiag[i];
447:     countA = ai[i+1] - adiag[i];
448:     countB = bi[i+1] - bi[i];
449:     bjj    = bj + bi[i];
450:     v2     = bv + bi[i];

452:      /* A-part */
453:     for (j=0; j<countA; j++){
454:       if (reuse == MAT_INITIAL_MATRIX) {
455:         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
456:       }
457:       val[jj++] = v1[j];
458:     }

460:     /* B-part */
461:     for(j=0; j < countB; j++){
462:       if (garray[bjj[j]] > rstart) {
463:         if (reuse == MAT_INITIAL_MATRIX) {
464:           row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
465:         }
466:         val[jj++] = v2[j];
467:       }
468:     }
469:     irow++;
470:   }
471:   return(0);
472: }

476: PetscErrorCode MatDestroy_MUMPS(Mat A)
477: {
478:   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;

482:   if (lu && lu->CleanUpMUMPS) {
483:     /* Terminate instance, deallocate memories */
484:     PetscFree2(lu->id.sol_loc,lu->id.isol_loc);
485:     VecScatterDestroy(&lu->scat_rhs);
486:     VecDestroy(&lu->b_seq);
487:     VecScatterDestroy(&lu->scat_sol);
488:     VecDestroy(&lu->x_seq);
489:     ierr=PetscFree(lu->id.perm_in);
490:     PetscFree(lu->irn);
491:     lu->id.job=JOB_END;
492: #if defined(PETSC_USE_COMPLEX)
493:     zmumps_c(&lu->id);
494: #else
495:     dmumps_c(&lu->id);
496: #endif
497:     MPI_Comm_free(&(lu->comm_mumps));
498:   }
499:   if (lu && lu->Destroy) {
500:     (lu->Destroy)(A);
501:   }
502:   PetscFree(A->spptr);

504:   /* clear composed functions */
505:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatFactorGetSolverPackage_C","",PETSC_NULL);
506:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMumpsSetIcntl_C","",PETSC_NULL);
507:   return(0);
508: }

512: PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x)
513: {
514:   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
515:   PetscScalar    *array;
516:   Vec            b_seq;
517:   IS             is_iden,is_petsc;
519:   PetscInt       i;

522:   lu->id.nrhs = 1;
523:   b_seq = lu->b_seq;
524:   if (lu->size > 1){
525:     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
526:     VecScatterBegin(lu->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);
527:     VecScatterEnd(lu->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);
528:     if (!lu->myid) {VecGetArray(b_seq,&array);}
529:   } else {  /* size == 1 */
530:     VecCopy(b,x);
531:     VecGetArray(x,&array);
532:   }
533:   if (!lu->myid) { /* define rhs on the host */
534:     lu->id.nrhs = 1;
535: #if defined(PETSC_USE_COMPLEX)
536:     lu->id.rhs = (mumps_double_complex*)array;
537: #else
538:     lu->id.rhs = array;
539: #endif
540:   }

542:   /* solve phase */
543:   /*-------------*/
544:   lu->id.job = JOB_SOLVE;
545: #if defined(PETSC_USE_COMPLEX)
546:   zmumps_c(&lu->id);
547: #else
548:   dmumps_c(&lu->id);
549: #endif
550:   if (lu->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1));

552:   if (lu->size > 1) { /* convert mumps distributed solution to petsc mpi x */
553:     if (!lu->nSolve){ /* create scatter scat_sol */
554:       ISCreateStride(PETSC_COMM_SELF,lu->id.lsol_loc,0,1,&is_iden); /* from */
555:       for (i=0; i<lu->id.lsol_loc; i++){
556:         lu->id.isol_loc[i] -= 1; /* change Fortran style to C style */
557:       }
558:       ISCreateGeneral(PETSC_COMM_SELF,lu->id.lsol_loc,lu->id.isol_loc,PETSC_COPY_VALUES,&is_petsc);  /* to */
559:       VecScatterCreate(lu->x_seq,is_iden,x,is_petsc,&lu->scat_sol);
560:       ISDestroy(&is_iden);
561:       ISDestroy(&is_petsc);
562:     }
563:     VecScatterBegin(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);
564:     VecScatterEnd(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);
565:   }
566:   lu->nSolve++;
567:   return(0);
568: }

572: PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x)
573: {
574:   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;

578:   lu->id.ICNTL(9) = 0;
579:   MatSolve_MUMPS(A,b,x);
580:   lu->id.ICNTL(9) = 1;
581:   return(0);
582: }

586: PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X)
587: {
589:   PetscBool      flg;

592:   PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,PETSC_NULL);
593:   if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
594:   PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,PETSC_NULL);
595:   if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");  SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_MUMPS() is not implemented yet");
596:   return(0);
597: }

599: #if !defined(PETSC_USE_COMPLEX)
600: /* 
601:   input:
602:    F:        numeric factor
603:   output:
604:    nneg:     total number of negative pivots
605:    nzero:    0
606:    npos:     (global dimension of F) - nneg
607: */

611: PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
612: {
613:   Mat_MUMPS      *lu =(Mat_MUMPS*)F->spptr;
615:   PetscMPIInt    size;

618:   MPI_Comm_size(((PetscObject)F)->comm,&size);
619:   /* 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 */
620:   if (size > 1 && lu->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",lu->id.INFOG(13));
621:   if (nneg){
622:     if (!lu->myid){
623:       *nneg = lu->id.INFOG(12);
624:     }
625:     MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);
626:   }
627:   if (nzero) *nzero = 0;
628:   if (npos)  *npos  = F->rmap->N - (*nneg);
629:   return(0);
630: }
631: #endif /* !defined(PETSC_USE_COMPLEX) */

635: PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
636: {
637:   Mat_MUMPS       *lu =(Mat_MUMPS*)(F)->spptr;
638:   PetscErrorCode  ierr;
639:   Mat             F_diag;
640:   PetscBool       isMPIAIJ;

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

645:   /* numerical factorization phase */
646:   /*-------------------------------*/
647:   lu->id.job = JOB_FACTNUMERIC;
648:   if(!lu->id.ICNTL(18)) {
649:     if (!lu->myid) {
650: #if defined(PETSC_USE_COMPLEX)
651:       lu->id.a = (mumps_double_complex*)lu->val;
652: #else
653:       lu->id.a = lu->val;
654: #endif
655:     }
656:   } else {
657: #if defined(PETSC_USE_COMPLEX)
658:     lu->id.a_loc = (mumps_double_complex*)lu->val;
659: #else
660:     lu->id.a_loc = lu->val;
661: #endif
662:   }
663: #if defined(PETSC_USE_COMPLEX)
664:   zmumps_c(&lu->id);
665: #else
666:   dmumps_c(&lu->id);
667: #endif
668:   if (lu->id.INFOG(1) < 0) {
669:     if (lu->id.INFO(1) == -13) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",lu->id.INFO(2));
670:     else SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: INFO(1)=%d, INFO(2)=%d\n",lu->id.INFO(1),lu->id.INFO(2));
671:   }
672:   if (!lu->myid && lu->id.ICNTL(16) > 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"  lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16));

674:   if (lu->size > 1){
675:     PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);
676:     if(isMPIAIJ) {
677:       F_diag = ((Mat_MPIAIJ *)(F)->data)->A;
678:     } else {
679:       F_diag = ((Mat_MPISBAIJ *)(F)->data)->A;
680:     }
681:     F_diag->assembled = PETSC_TRUE;
682:     if (lu->nSolve){
683:       VecScatterDestroy(&lu->scat_sol);
684:       PetscFree2(lu->id.sol_loc,lu->id.isol_loc);
685:       VecDestroy(&lu->x_seq);
686:     }
687:   }
688:   (F)->assembled   = PETSC_TRUE;
689:   lu->matstruc     = SAME_NONZERO_PATTERN;
690:   lu->CleanUpMUMPS = PETSC_TRUE;
691:   lu->nSolve       = 0;
692: 
693:   if (lu->size > 1){
694:     /* distributed solution */
695:     if (!lu->nSolve){
696:       /* Create x_seq=sol_loc for repeated use */
697:       PetscInt    lsol_loc;
698:       PetscScalar *sol_loc;
699:       lsol_loc = lu->id.INFO(23); /* length of sol_loc */
700:       PetscMalloc2(lsol_loc,PetscScalar,&sol_loc,lsol_loc,PetscInt,&lu->id.isol_loc);
701:       lu->id.lsol_loc = lsol_loc;
702: #if defined(PETSC_USE_COMPLEX)
703:       lu->id.sol_loc  = (mumps_double_complex*)sol_loc;
704: #else
705:       lu->id.sol_loc  = sol_loc;
706: #endif
707:       VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&lu->x_seq);
708:     }
709:   }
710:   return(0);
711: }

713: /* Sets MUMPS options from the options database */
716: PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A)
717: {
718:   Mat_MUMPS        *mumps = (Mat_MUMPS*)F->spptr;
719:   PetscErrorCode   ierr;
720:   PetscInt         icntl;
721:   PetscBool        flg;

724:   PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"MUMPS Options","Mat");
725:   PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);
726:   if (flg) mumps->id.ICNTL(1) = icntl;
727:   PetscOptionsInt("-mat_mumps_icntl_2","ICNTL(2): output stream for diagnostic printing, statistics, and warning","None",mumps->id.ICNTL(2),&icntl,&flg);
728:   if (flg) mumps->id.ICNTL(2) = icntl;
729:   PetscOptionsInt("-mat_mumps_icntl_3","ICNTL(3): output stream for global information, collected on the host","None",mumps->id.ICNTL(3),&icntl,&flg);
730:   if (flg) mumps->id.ICNTL(3) = icntl;

732:   PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",mumps->id.ICNTL(4),&icntl,&flg);
733:   if (flg) mumps->id.ICNTL(4) = icntl;
734:   if (mumps->id.ICNTL(4) || PetscLogPrintInfo ) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */
735: 
736:   PetscOptionsInt("-mat_mumps_icntl_6","ICNTL(6): permuting and/or scaling the matrix (0 to 7)","None",mumps->id.ICNTL(6),&icntl,&flg);
737:   if (flg) mumps->id.ICNTL(6) = icntl;

739:   PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7). 3=Scotch, 4=PORD, 5=Metis","None",mumps->id.ICNTL(7),&icntl,&flg);
740:   if (flg) {
741:     if (icntl== 1 && mumps->size > 1){
742:       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");
743:     } else {
744:       mumps->id.ICNTL(7) = icntl;
745:     }
746:   }
747: 
748:   PetscOptionsInt("-mat_mumps_icntl_8","ICNTL(8): scaling strategy (-2 to 8 or 77)","None",mumps->id.ICNTL(8),&mumps->id.ICNTL(8),PETSC_NULL);
749:   PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",mumps->id.ICNTL(10),&mumps->id.ICNTL(10),PETSC_NULL);
750:   PetscOptionsInt("-mat_mumps_icntl_11","ICNTL(11): statistics related to the linear system solved (via -ksp_view)","None",mumps->id.ICNTL(11),&mumps->id.ICNTL(11),PETSC_NULL);
751:   PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control: defines the ordering strategy with scaling constraints (0 to 3)","None",mumps->id.ICNTL(12),&mumps->id.ICNTL(12),PETSC_NULL);
752:   PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control: with or without ScaLAPACK","None",mumps->id.ICNTL(13),&mumps->id.ICNTL(13),PETSC_NULL);
753:   PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage of estimated workspace increase","None",mumps->id.ICNTL(14),&mumps->id.ICNTL(14),PETSC_NULL);
754:   PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): Schur complement","None",mumps->id.ICNTL(19),&mumps->id.ICNTL(19),PETSC_NULL);

756:   PetscOptionsInt("-mat_mumps_icntl_22","ICNTL(22): in-core/out-of-core facility (0 or 1)","None",mumps->id.ICNTL(22),&mumps->id.ICNTL(22),PETSC_NULL);
757:   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),PETSC_NULL);
758:   PetscOptionsInt("-mat_mumps_icntl_24","ICNTL(24): detection of null pivot rows (0 or 1)","None",mumps->id.ICNTL(24),&mumps->id.ICNTL(24),PETSC_NULL);
759:   if (mumps->id.ICNTL(24)){
760:     mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */
761:   }

763:   PetscOptionsInt("-mat_mumps_icntl_25","ICNTL(25): computation of a null space basis","None",mumps->id.ICNTL(25),&mumps->id.ICNTL(25),PETSC_NULL);
764:   PetscOptionsInt("-mat_mumps_icntl_26","ICNTL(26): Schur options for right-hand side or solution vector","None",mumps->id.ICNTL(26),&mumps->id.ICNTL(26),PETSC_NULL);
765:   PetscOptionsInt("-mat_mumps_icntl_27","ICNTL(27): experimental parameter","None",mumps->id.ICNTL(27),&mumps->id.ICNTL(27),PETSC_NULL);
766:   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),PETSC_NULL);
767:   PetscOptionsInt("-mat_mumps_icntl_29","ICNTL(29): parallel ordering 1 = ptscotch 2 = parmetis","None",mumps->id.ICNTL(29),&mumps->id.ICNTL(29),PETSC_NULL);
768:   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),PETSC_NULL);
769:   PetscOptionsInt("-mat_mumps_icntl_31","ICNTL(31): factors can be discarded in the solve phase","None",mumps->id.ICNTL(31),&mumps->id.ICNTL(31),PETSC_NULL);
770:   PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),PETSC_NULL);

772:   PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),PETSC_NULL);
773:   PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",mumps->id.CNTL(2),&mumps->id.CNTL(2),PETSC_NULL);
774:   PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),PETSC_NULL);
775:   PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",mumps->id.CNTL(4),&mumps->id.CNTL(4),PETSC_NULL);
776:   PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",mumps->id.CNTL(5),&mumps->id.CNTL(5),PETSC_NULL);
777: 
778:   PetscOptionsString("-mat_mumps_ooc_tmpdir", "out of core directory", "None", mumps->id.ooc_tmpdir, mumps->id.ooc_tmpdir, 256, PETSC_NULL);
779:   PetscOptionsEnd();
780:   return(0);
781: }
782: 
785: PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS* mumps)
786: {
787:   PetscErrorCode  ierr;

790:   MPI_Comm_rank(((PetscObject)A)->comm, &mumps->myid);
791:   MPI_Comm_size(((PetscObject)A)->comm,&mumps->size);
792:   MPI_Comm_dup(((PetscObject)A)->comm,&(mumps->comm_mumps));
793:   mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps);

795:   mumps->id.job = JOB_INIT;
796:   mumps->id.par = 1;  /* host participates factorizaton and solve */
797:   mumps->id.sym = mumps->sym;
798: #if defined(PETSC_USE_COMPLEX)
799:   zmumps_c(&mumps->id);
800: #else
801:   dmumps_c(&mumps->id);
802: #endif

804:   mumps->CleanUpMUMPS = PETSC_FALSE;
805:   mumps->scat_rhs     = PETSC_NULL;
806:   mumps->scat_sol     = PETSC_NULL;
807:   mumps->nSolve       = 0;

809:   /* set PETSc-MUMPS default options - override MUMPS default */
810:   mumps->id.ICNTL(3) = 0;
811:   mumps->id.ICNTL(4) = 0;
812:   if (mumps->size == 1){
813:     mumps->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
814:   } else {
815:     mumps->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
816:     mumps->id.ICNTL(21) = 1;   /* distributed solution */
817:   }
818:   return(0);
819: }
820: 
821: /* Note the Petsc r and c permutations are ignored */
824: PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
825: {
826:   Mat_MUMPS          *lu = (Mat_MUMPS*)F->spptr;
827:   PetscErrorCode     ierr;
828:   Vec                b;
829:   IS                 is_iden;
830:   const PetscInt     M = A->rmap->N;

833:   lu->matstruc = DIFFERENT_NONZERO_PATTERN;

835:   /* Set MUMPS options from the options database */
836:   PetscSetMUMPSFromOptions(F,A);
837: 
838:   (*lu->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &lu->nz, &lu->irn, &lu->jcn, &lu->val);

840:   /* analysis phase */
841:   /*----------------*/
842:   lu->id.job = JOB_FACTSYMBOLIC;
843:   lu->id.n = M;
844:   switch (lu->id.ICNTL(18)){
845:   case 0:  /* centralized assembled matrix input */
846:     if (!lu->myid) {
847:       lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
848:       if (lu->id.ICNTL(6)>1){
849: #if defined(PETSC_USE_COMPLEX)
850:         lu->id.a = (mumps_double_complex*)lu->val;
851: #else
852:         lu->id.a = lu->val;
853: #endif
854:       }
855:       if (lu->id.ICNTL(7) == 1){ /* use user-provide matrix ordering */
856:         if (!lu->myid) {
857:           const PetscInt *idx;
858:           PetscInt i,*perm_in;
859:           PetscMalloc(M*sizeof(PetscInt),&perm_in);
860:           ISGetIndices(r,&idx);
861:           lu->id.perm_in = perm_in;
862:           for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */
863:           ISRestoreIndices(r,&idx);
864:         }
865:       }
866:     }
867:     break;
868:   case 3:  /* distributed assembled matrix input (size>1) */
869:     lu->id.nz_loc = lu->nz;
870:     lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
871:     if (lu->id.ICNTL(6)>1) {
872: #if defined(PETSC_USE_COMPLEX)
873:       lu->id.a_loc = (mumps_double_complex*)lu->val;
874: #else
875:       lu->id.a_loc = lu->val;
876: #endif
877:     }
878:     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
879:     if (!lu->myid){
880:       VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);
881:       ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);
882:     } else {
883:       VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);
884:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);
885:     }
886:     VecCreate(((PetscObject)A)->comm,&b);
887:     VecSetSizes(b,A->rmap->n,PETSC_DECIDE);
888:     VecSetFromOptions(b);

890:     VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);
891:     ISDestroy(&is_iden);
892:     VecDestroy(&b);
893:     break;
894:     }
895: #if defined(PETSC_USE_COMPLEX)
896:   zmumps_c(&lu->id);
897: #else
898:   dmumps_c(&lu->id);
899: #endif
900:   if (lu->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
901: 
902:   F->ops->lufactornumeric  = MatFactorNumeric_MUMPS;
903:   F->ops->solve            = MatSolve_MUMPS;
904:   F->ops->solvetranspose   = MatSolveTranspose_MUMPS;
905:   F->ops->matsolve         = MatMatSolve_MUMPS;
906:   return(0);
907: }

909: /* Note the Petsc r and c permutations are ignored */
912: PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
913: {

915:   Mat_MUMPS       *lu = (Mat_MUMPS*)F->spptr;
916:   PetscErrorCode  ierr;
917:   Vec             b;
918:   IS              is_iden;
919:   const PetscInt  M = A->rmap->N;

922:   lu->matstruc = DIFFERENT_NONZERO_PATTERN;

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

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

929:   /* analysis phase */
930:   /*----------------*/
931:   lu->id.job = JOB_FACTSYMBOLIC;
932:   lu->id.n = M;
933:   switch (lu->id.ICNTL(18)){
934:   case 0:  /* centralized assembled matrix input */
935:     if (!lu->myid) {
936:       lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
937:       if (lu->id.ICNTL(6)>1){
938: #if defined(PETSC_USE_COMPLEX)
939:         lu->id.a = (mumps_double_complex*)lu->val;
940: #else
941:         lu->id.a = lu->val;
942: #endif
943:       }
944:     }
945:     break;
946:   case 3:  /* distributed assembled matrix input (size>1) */
947:     lu->id.nz_loc = lu->nz;
948:     lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
949:     if (lu->id.ICNTL(6)>1) {
950: #if defined(PETSC_USE_COMPLEX)
951:       lu->id.a_loc = (mumps_double_complex*)lu->val;
952: #else
953:       lu->id.a_loc = lu->val;
954: #endif
955:     }
956:     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
957:     if (!lu->myid){
958:       VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);
959:       ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);
960:     } else {
961:       VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);
962:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);
963:     }
964:     VecCreate(((PetscObject)A)->comm,&b);
965:     VecSetSizes(b,A->rmap->n,PETSC_DECIDE);
966:     VecSetFromOptions(b);

968:     VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);
969:     ISDestroy(&is_iden);
970:     VecDestroy(&b);
971:     break;
972:     }
973: #if defined(PETSC_USE_COMPLEX)
974:   zmumps_c(&lu->id);
975: #else
976:   dmumps_c(&lu->id);
977: #endif
978:   if (lu->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
979: 
980:   F->ops->lufactornumeric  = MatFactorNumeric_MUMPS;
981:   F->ops->solve            = MatSolve_MUMPS;
982:   F->ops->solvetranspose   = MatSolveTranspose_MUMPS;
983:   return(0);
984: }

986: /* Note the Petsc r permutation and factor info are ignored */
989: PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
990: {
991:   Mat_MUMPS          *lu = (Mat_MUMPS*)F->spptr;
992:   PetscErrorCode     ierr;
993:   Vec                b;
994:   IS                 is_iden;
995:   const PetscInt     M = A->rmap->N;

998:   lu->matstruc = DIFFERENT_NONZERO_PATTERN;

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

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

1005:   /* analysis phase */
1006:   /*----------------*/
1007:   lu->id.job = JOB_FACTSYMBOLIC;
1008:   lu->id.n = M;
1009:   switch (lu->id.ICNTL(18)){
1010:   case 0:  /* centralized assembled matrix input */
1011:     if (!lu->myid) {
1012:       lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
1013:       if (lu->id.ICNTL(6)>1){
1014: #if defined(PETSC_USE_COMPLEX)
1015:         lu->id.a = (mumps_double_complex*)lu->val;
1016: #else
1017:         lu->id.a = lu->val;
1018: #endif
1019:       }
1020:     }
1021:     break;
1022:   case 3:  /* distributed assembled matrix input (size>1) */
1023:     lu->id.nz_loc = lu->nz;
1024:     lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
1025:     if (lu->id.ICNTL(6)>1) {
1026: #if defined(PETSC_USE_COMPLEX)
1027:       lu->id.a_loc = (mumps_double_complex*)lu->val;
1028: #else
1029:       lu->id.a_loc = lu->val;
1030: #endif
1031:     }
1032:     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1033:     if (!lu->myid){
1034:       VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);
1035:       ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);
1036:     } else {
1037:       VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);
1038:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);
1039:     }
1040:     VecCreate(((PetscObject)A)->comm,&b);
1041:     VecSetSizes(b,A->rmap->n,PETSC_DECIDE);
1042:     VecSetFromOptions(b);

1044:     VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);
1045:     ISDestroy(&is_iden);
1046:     VecDestroy(&b);
1047:     break;
1048:     }
1049: #if defined(PETSC_USE_COMPLEX)
1050:   zmumps_c(&lu->id);
1051: #else
1052:   dmumps_c(&lu->id);
1053: #endif
1054:   if (lu->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));

1056:   F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
1057:   F->ops->solve                 = MatSolve_MUMPS;
1058:   F->ops->solvetranspose        = MatSolve_MUMPS;
1059:   F->ops->matsolve              = MatMatSolve_MUMPS;
1060: #if !defined(PETSC_USE_COMPLEX)
1061:   F->ops->getinertia            = MatGetInertia_SBAIJMUMPS;
1062: #else
1063:   F->ops->getinertia            = PETSC_NULL;
1064: #endif
1065:   return(0);
1066: }

1070: PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
1071: {
1072:   PetscErrorCode    ierr;
1073:   PetscBool         iascii;
1074:   PetscViewerFormat format;
1075:   Mat_MUMPS         *lu=(Mat_MUMPS*)A->spptr;

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

1081:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1082:   if (iascii) {
1083:     PetscViewerGetFormat(viewer,&format);
1084:     if (format == PETSC_VIEWER_ASCII_INFO){
1085:       PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");
1086:       PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                   %d \n",lu->id.sym);
1087:       PetscViewerASCIIPrintf(viewer,"  PAR (host participation):            %d \n",lu->id.par);
1088:       PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):         %d \n",lu->id.ICNTL(1));
1089:       PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg): %d \n",lu->id.ICNTL(2));
1090:       PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):   %d \n",lu->id.ICNTL(3));
1091:       PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):        %d \n",lu->id.ICNTL(4));
1092:       PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):         %d \n",lu->id.ICNTL(5));
1093:       PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):        %d \n",lu->id.ICNTL(6));
1094:       PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (sequentia matrix ordering):%d \n",lu->id.ICNTL(7));
1095:       PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scalling strategy):        %d \n",lu->id.ICNTL(8));
1096:       PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements):  %d \n",lu->id.ICNTL(10));
1097:       PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):          %d \n",lu->id.ICNTL(11));
1098:       if (lu->id.ICNTL(11)>0) {
1099:         PetscViewerASCIIPrintf(viewer,"    RINFOG(4) (inf norm of input mat):        %g\n",lu->id.RINFOG(4));
1100:         PetscViewerASCIIPrintf(viewer,"    RINFOG(5) (inf norm of solution):         %g\n",lu->id.RINFOG(5));
1101:         PetscViewerASCIIPrintf(viewer,"    RINFOG(6) (inf norm of residual):         %g\n",lu->id.RINFOG(6));
1102:         PetscViewerASCIIPrintf(viewer,"    RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",lu->id.RINFOG(7),lu->id.RINFOG(8));
1103:         PetscViewerASCIIPrintf(viewer,"    RINFOG(9) (error estimate):               %g \n",lu->id.RINFOG(9));
1104:         PetscViewerASCIIPrintf(viewer,"    RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",lu->id.RINFOG(10),lu->id.RINFOG(11));
1105:       }
1106:       PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",lu->id.ICNTL(12));
1107:       PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",lu->id.ICNTL(13));
1108:       PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));
1109:       /* ICNTL(15-17) not used */
1110:       PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",lu->id.ICNTL(18));
1111:       PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Shur complement info):                       %d \n",lu->id.ICNTL(19));
1112:       PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",lu->id.ICNTL(20));
1113:       PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d \n",lu->id.ICNTL(21));
1114:       PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",lu->id.ICNTL(22));
1115:       PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",lu->id.ICNTL(23));
1116: 
1117:       PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",lu->id.ICNTL(24));
1118:       PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",lu->id.ICNTL(25));
1119:       PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",lu->id.ICNTL(26));
1120:       PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",lu->id.ICNTL(27));
1121:       PetscViewerASCIIPrintf(viewer,"  ICNTL(28) (use parallel or sequential ordering):        %d \n",lu->id.ICNTL(28));
1122:       PetscViewerASCIIPrintf(viewer,"  ICNTL(29) (parallel ordering):                          %d \n",lu->id.ICNTL(29));
1123: 
1124:       PetscViewerASCIIPrintf(viewer,"  ICNTL(30) (user-specified set of entries in inv(A)):    %d \n",lu->id.ICNTL(30));
1125:       PetscViewerASCIIPrintf(viewer,"  ICNTL(31) (factors is discarded in the solve phase):    %d \n",lu->id.ICNTL(31));
1126:       PetscViewerASCIIPrintf(viewer,"  ICNTL(33) (compute determinant):                        %d \n",lu->id.ICNTL(33));
1127: 
1128:       PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",lu->id.CNTL(1));
1129:       PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));
1130:       PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",lu->id.CNTL(3));
1131:       PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",lu->id.CNTL(4));
1132:       PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",lu->id.CNTL(5));
1133: 
1134:       /* infomation local to each processor */
1135:       PetscViewerASCIIPrintf(viewer, "  RINFO(1) (local estimated flops for the elimination after analysis): \n");
1136:       PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);
1137:       PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %g \n",lu->myid,lu->id.RINFO(1));
1138:       PetscViewerFlush(viewer);
1139:       PetscViewerASCIIPrintf(viewer, "  RINFO(2) (local estimated flops for the assembly after factorization): \n");
1140:       PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",lu->myid,lu->id.RINFO(2));
1141:       PetscViewerFlush(viewer);
1142:       PetscViewerASCIIPrintf(viewer, "  RINFO(3) (local estimated flops for the elimination after factorization): \n");
1143:       PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",lu->myid,lu->id.RINFO(3));
1144:       PetscViewerFlush(viewer);
1145: 
1146:       PetscViewerASCIIPrintf(viewer, "  INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");
1147:       PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] %d \n",lu->myid,lu->id.INFO(15));
1148:       PetscViewerFlush(viewer);
1149: 
1150:       PetscViewerASCIIPrintf(viewer, "  INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");
1151:       PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",lu->myid,lu->id.INFO(16));
1152:       PetscViewerFlush(viewer);
1153: 
1154:       PetscViewerASCIIPrintf(viewer, "  INFO(23) (num of pivots eliminated on this processor after factorization): \n");
1155:       PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",lu->myid,lu->id.INFO(23));
1156:       PetscViewerFlush(viewer);
1157:       PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);

1159:       if (!lu->myid){ /* information from the host */
1160:         PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));
1161:         PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));
1162:         PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));
1163:         PetscViewerASCIIPrintf(viewer,"  (RINFOG(12) RINFOG(13))*2^INFOG(34) (determinant): (%g,%g)*(2^%d)\n",lu->id.RINFOG(12),lu->id.RINFOG(13),lu->id.INFOG(34));
1164: 
1165:         PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));
1166:         PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));
1167:         PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));
1168:         PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));
1169:         PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively use after analysis): %d \n",lu->id.INFOG(7));
1170:         PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));
1171:         PetscViewerASCIIPrintf(viewer,"  INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",lu->id.INFOG(9));
1172:         PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",lu->id.INFOG(10));
1173:         PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",lu->id.INFOG(11));
1174:         PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));
1175:         PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));
1176:         PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));
1177:         PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));
1178:         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",lu->id.INFOG(16));
1179:         PetscViewerASCIIPrintf(viewer,"  INFOG(17) (estimated size of all MUMPS internal data for factorization after analysis: sum over all processors): %d \n",lu->id.INFOG(17));
1180:         PetscViewerASCIIPrintf(viewer,"  INFOG(18) (size of all MUMPS internal data allocated during factorization: value on the most memory consuming processor): %d \n",lu->id.INFOG(18));
1181:         PetscViewerASCIIPrintf(viewer,"  INFOG(19) (size of all MUMPS internal data allocated during factorization: sum over all processors): %d \n",lu->id.INFOG(19));
1182:         PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));
1183:         PetscViewerASCIIPrintf(viewer,"  INFOG(21) (size in MB of memory effectively used during factorization - value on the most memory consuming processor): %d \n",lu->id.INFOG(21));
1184:         PetscViewerASCIIPrintf(viewer,"  INFOG(22) (size in MB of memory effectively used during factorization - sum over all processors): %d \n",lu->id.INFOG(22));
1185:         PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",lu->id.INFOG(23));
1186:         PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",lu->id.INFOG(24));
1187:         PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",lu->id.INFOG(25));
1188:       }
1189:     }
1190:   }
1191:   return(0);
1192: }

1196: PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
1197: {
1198:   Mat_MUMPS  *mumps =(Mat_MUMPS*)A->spptr;

1201:   info->block_size        = 1.0;
1202:   info->nz_allocated      = mumps->id.INFOG(20);
1203:   info->nz_used           = mumps->id.INFOG(20);
1204:   info->nz_unneeded       = 0.0;
1205:   info->assemblies        = 0.0;
1206:   info->mallocs           = 0.0;
1207:   info->memory            = 0.0;
1208:   info->fill_ratio_given  = 0;
1209:   info->fill_ratio_needed = 0;
1210:   info->factor_mallocs    = 0;
1211:   return(0);
1212: }

1214: /* -------------------------------------------------------------------------------------------*/
1217: PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival)
1218: {
1219:   Mat_MUMPS *lu =(Mat_MUMPS*)F->spptr;

1222:   lu->id.ICNTL(icntl) = ival;
1223:   return(0);
1224: }

1228: /*@
1229:   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()

1231:    Logically Collective on Mat

1233:    Input Parameters:
1234: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1235: .  icntl - index of MUMPS parameter array ICNTL()
1236: -  ival - value of MUMPS ICNTL(icntl)

1238:   Options Database:
1239: .   -mat_mumps_icntl_<icntl> <ival>

1241:    Level: beginner

1243:    References: MUMPS Users' Guide 

1245: .seealso: MatGetFactor()
1246: @*/
1247: PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival)
1248: {

1254:   PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));
1255:   return(0);
1256: }

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

1262:   Works with MATAIJ and MATSBAIJ matrices

1264:   Options Database Keys:
1265: + -mat_mumps_icntl_4 <0,...,4> - print level
1266: . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
1267: . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guidec)
1268: . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
1269: . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
1270: . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
1271: . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
1272: . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
1273: . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
1274: . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
1275: . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
1276: . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
1277: - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold

1279:   Level: beginner

1281: .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage

1283: M*/

1285: EXTERN_C_BEGIN
1288: PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type)
1289: {
1291:   *type = MATSOLVERMUMPS;
1292:   return(0);
1293: }
1294: EXTERN_C_END

1296: EXTERN_C_BEGIN
1297: /* MatGetFactor for Seq and MPI AIJ matrices */
1300: PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
1301: {
1302:   Mat            B;
1304:   Mat_MUMPS      *mumps;
1305:   PetscBool      isSeqAIJ;

1308:   /* Create the factorization matrix */
1309:   PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);
1310:   MatCreate(((PetscObject)A)->comm,&B);
1311:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
1312:   MatSetType(B,((PetscObject)A)->type_name);
1313:   if (isSeqAIJ) {
1314:     MatSeqAIJSetPreallocation(B,0,PETSC_NULL);
1315:   } else {
1316:     MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);
1317:   }

1319:   PetscNewLog(B,Mat_MUMPS,&mumps);
1320:   B->ops->view             = MatView_MUMPS;
1321:   B->ops->getinfo          = MatGetInfo_MUMPS;
1322:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);
1323:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl_MUMPS",MatMumpsSetIcntl_MUMPS);
1324:   if (ftype == MAT_FACTOR_LU) {
1325:     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
1326:     B->factortype = MAT_FACTOR_LU;
1327:     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
1328:     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
1329:     mumps->sym = 0;
1330:   } else {
1331:     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
1332:     B->factortype = MAT_FACTOR_CHOLESKY;
1333:     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
1334:     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
1335:     if (A->spd_set && A->spd) mumps->sym = 1;
1336:     else                      mumps->sym = 2;
1337:   }

1339:   mumps->isAIJ        = PETSC_TRUE;
1340:   mumps->Destroy      = B->ops->destroy;
1341:   B->ops->destroy     = MatDestroy_MUMPS;
1342:   B->spptr            = (void*)mumps;
1343:   PetscInitializeMUMPS(A,mumps);

1345:   *F = B;
1346:   return(0);
1347: }
1348: EXTERN_C_END


1351: EXTERN_C_BEGIN
1352: /* MatGetFactor for Seq and MPI SBAIJ matrices */
1355: PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
1356: {
1357:   Mat            B;
1359:   Mat_MUMPS      *mumps;
1360:   PetscBool      isSeqSBAIJ;

1363:   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
1364:   if(A->rmap->bs > 1) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with block size > 1 with MUMPS Cholesky, use AIJ matrix instead");
1365:   PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);
1366:   /* Create the factorization matrix */
1367:   MatCreate(((PetscObject)A)->comm,&B);
1368:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
1369:   MatSetType(B,((PetscObject)A)->type_name);
1370:   PetscNewLog(B,Mat_MUMPS,&mumps);
1371:   if (isSeqSBAIJ) {
1372:     MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);
1373:     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
1374:   } else {
1375:     MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);
1376:     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
1377:   }

1379:   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
1380:   B->ops->view                   = MatView_MUMPS;
1381:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);
1382:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl",MatMumpsSetIcntl);
1383:   B->factortype                  = MAT_FACTOR_CHOLESKY;
1384:   if (A->spd_set && A->spd) mumps->sym = 1;
1385:   else                      mumps->sym = 2;

1387:   mumps->isAIJ        = PETSC_FALSE;
1388:   mumps->Destroy      = B->ops->destroy;
1389:   B->ops->destroy     = MatDestroy_MUMPS;
1390:   B->spptr            = (void*)mumps;
1391:   PetscInitializeMUMPS(A,mumps);

1393:   *F = B;
1394:   return(0);
1395: }
1396: EXTERN_C_END

1398: EXTERN_C_BEGIN
1401: PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
1402: {
1403:   Mat            B;
1405:   Mat_MUMPS      *mumps;
1406:   PetscBool      isSeqBAIJ;

1409:   /* Create the factorization matrix */
1410:   PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);
1411:   MatCreate(((PetscObject)A)->comm,&B);
1412:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
1413:   MatSetType(B,((PetscObject)A)->type_name);
1414:   if (isSeqBAIJ) {
1415:     MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,PETSC_NULL);
1416:   } else {
1417:     MatMPIBAIJSetPreallocation(B,A->rmap->bs,0,PETSC_NULL,0,PETSC_NULL);
1418:   }

1420:   PetscNewLog(B,Mat_MUMPS,&mumps);
1421:   if (ftype == MAT_FACTOR_LU) {
1422:     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
1423:     B->factortype = MAT_FACTOR_LU;
1424:     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
1425:     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
1426:     mumps->sym = 0;
1427:   } else {
1428:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n");
1429:   }

1431:   B->ops->view             = MatView_MUMPS;
1432:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);
1433:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl_MUMPS",MatMumpsSetIcntl_MUMPS);

1435:   mumps->isAIJ        = PETSC_TRUE;
1436:   mumps->Destroy      = B->ops->destroy;
1437:   B->ops->destroy     = MatDestroy_MUMPS;
1438:   B->spptr            = (void*)mumps;
1439:   PetscInitializeMUMPS(A,mumps);

1441:   *F = B;
1442:   return(0);
1443: }
1444: EXTERN_C_END