xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision 38be02ac5e7f436cf326ff72e1b1b2b59afbca4d)
11c2a3de1SBarry Smith 
2397b6df1SKris Buschelman /*
3c2b5dc30SHong Zhang     Provides an interface to the MUMPS sparse solver
4397b6df1SKris Buschelman */
551d5961aSHong Zhang 
6c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> /*I  "petscmat.h"  I*/
7c6db04a5SJed Brown #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
8b5fa320bSStefano Zampini #include <petscblaslapack.h>
9397b6df1SKris Buschelman 
10397b6df1SKris Buschelman EXTERN_C_BEGIN
11397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
122907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
132907cef9SHong Zhang #include <cmumps_c.h>
142907cef9SHong Zhang #else
15c6db04a5SJed Brown #include <zmumps_c.h>
162907cef9SHong Zhang #endif
172907cef9SHong Zhang #else
182907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
192907cef9SHong Zhang #include <smumps_c.h>
20397b6df1SKris Buschelman #else
21c6db04a5SJed Brown #include <dmumps_c.h>
22397b6df1SKris Buschelman #endif
232907cef9SHong Zhang #endif
24397b6df1SKris Buschelman EXTERN_C_END
25397b6df1SKris Buschelman #define JOB_INIT -1
263d472b54SHong Zhang #define JOB_FACTSYMBOLIC 1
273d472b54SHong Zhang #define JOB_FACTNUMERIC 2
283d472b54SHong Zhang #define JOB_SOLVE 3
29397b6df1SKris Buschelman #define JOB_END -2
303d472b54SHong Zhang 
312907cef9SHong Zhang /* calls to MUMPS */
322907cef9SHong Zhang #if defined(PETSC_USE_COMPLEX)
332907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
342907cef9SHong Zhang #define PetscMUMPS_c cmumps_c
352907cef9SHong Zhang #else
362907cef9SHong Zhang #define PetscMUMPS_c zmumps_c
372907cef9SHong Zhang #endif
382907cef9SHong Zhang #else
392907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
402907cef9SHong Zhang #define PetscMUMPS_c smumps_c
412907cef9SHong Zhang #else
422907cef9SHong Zhang #define PetscMUMPS_c dmumps_c
432907cef9SHong Zhang #endif
442907cef9SHong Zhang #endif
452907cef9SHong Zhang 
46940cd9d6SSatish Balay /* declare MumpsScalar */
47940cd9d6SSatish Balay #if defined(PETSC_USE_COMPLEX)
48940cd9d6SSatish Balay #if defined(PETSC_USE_REAL_SINGLE)
49940cd9d6SSatish Balay #define MumpsScalar mumps_complex
50940cd9d6SSatish Balay #else
51940cd9d6SSatish Balay #define MumpsScalar mumps_double_complex
52940cd9d6SSatish Balay #endif
53940cd9d6SSatish Balay #else
54940cd9d6SSatish Balay #define MumpsScalar PetscScalar
55940cd9d6SSatish Balay #endif
563d472b54SHong Zhang 
57397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */
58397b6df1SKris Buschelman #define ICNTL(I) icntl[(I)-1]
59397b6df1SKris Buschelman #define CNTL(I) cntl[(I)-1]
60397b6df1SKris Buschelman #define INFOG(I) infog[(I)-1]
61a7aca84bSHong Zhang #define INFO(I) info[(I)-1]
62397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1]
63adc1d99fSHong Zhang #define RINFO(I) rinfo[(I)-1]
64397b6df1SKris Buschelman 
65397b6df1SKris Buschelman typedef struct {
66397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
672907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
682907cef9SHong Zhang   CMUMPS_STRUC_C id;
692907cef9SHong Zhang #else
70397b6df1SKris Buschelman   ZMUMPS_STRUC_C id;
712907cef9SHong Zhang #endif
722907cef9SHong Zhang #else
732907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
742907cef9SHong Zhang   SMUMPS_STRUC_C id;
75397b6df1SKris Buschelman #else
76397b6df1SKris Buschelman   DMUMPS_STRUC_C id;
77397b6df1SKris Buschelman #endif
782907cef9SHong Zhang #endif
792907cef9SHong Zhang 
80397b6df1SKris Buschelman   MatStructure matstruc;
81c1490034SHong Zhang   PetscMPIInt  myid,size;
82a5e57a09SHong Zhang   PetscInt     *irn,*jcn,nz,sym;
83397b6df1SKris Buschelman   PetscScalar  *val;
84397b6df1SKris Buschelman   MPI_Comm     comm_mumps;
856f3cc6f9SBarry Smith   PetscBool    isAIJ;
86a5e57a09SHong Zhang   PetscInt     ICNTL9_pre;           /* check if ICNTL(9) is changed from previous MatSolve */
87801fbe65SHong Zhang   VecScatter   scat_rhs, scat_sol;   /* used by MatSolve() */
88801fbe65SHong Zhang   Vec          b_seq,x_seq;
89b34f08ffSHong Zhang   PetscInt     ninfo,*info;          /* display INFO */
90b5fa320bSStefano Zampini   PetscInt     sizeredrhs;
91b5fa320bSStefano Zampini   PetscInt     *schur_pivots;
9259ac8732SStefano Zampini   PetscInt     schur_B_lwork;
93b5fa320bSStefano Zampini   PetscScalar  *schur_work;
9459ac8732SStefano Zampini   PetscScalar  *schur_sol;
9559ac8732SStefano Zampini   PetscInt     schur_sizesol;
9659ac8732SStefano Zampini   PetscBool    schur_factored;
9759ac8732SStefano Zampini   PetscBool    schur_inverted;
98a0b0af32SStefano Zampini   PetscInt     schur_sym;
992205254eSKarl Rupp 
100bccb9932SShri Abhyankar   PetscErrorCode (*ConvertToTriples)(Mat, int, MatReuse, int*, int**, int**, PetscScalar**);
101f0c56d0fSKris Buschelman } Mat_MUMPS;
102f0c56d0fSKris Buschelman 
10309573ac7SBarry Smith extern PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);
104b24902e0SBarry Smith 
10559ac8732SStefano Zampini static PetscErrorCode MatMumpsResetSchur_Private(Mat_MUMPS* mumps)
106b5fa320bSStefano Zampini {
107b5fa320bSStefano Zampini   PetscErrorCode ierr;
108b5fa320bSStefano Zampini 
109b5fa320bSStefano Zampini   PetscFunctionBegin;
11059ac8732SStefano Zampini   ierr = PetscFree2(mumps->id.listvar_schur,mumps->id.schur);CHKERRQ(ierr);
11159ac8732SStefano Zampini   ierr = PetscFree(mumps->id.redrhs);CHKERRQ(ierr);
11259ac8732SStefano Zampini   ierr = PetscFree(mumps->schur_sol);CHKERRQ(ierr);
11359ac8732SStefano Zampini   ierr = PetscFree(mumps->schur_pivots);CHKERRQ(ierr);
11459ac8732SStefano Zampini   ierr = PetscFree(mumps->schur_work);CHKERRQ(ierr);
11559ac8732SStefano Zampini   mumps->id.size_schur = 0;
11659ac8732SStefano Zampini   mumps->id.ICNTL(19) = 0;
11759ac8732SStefano Zampini   PetscFunctionReturn(0);
11859ac8732SStefano Zampini }
11959ac8732SStefano Zampini 
12059ac8732SStefano Zampini static PetscErrorCode MatMumpsFactorSchur_Private(Mat_MUMPS* mumps)
12159ac8732SStefano Zampini {
12259ac8732SStefano Zampini   PetscBLASInt   B_N,B_ierr,B_slda;
12359ac8732SStefano Zampini   PetscErrorCode ierr;
12459ac8732SStefano Zampini 
12559ac8732SStefano Zampini   PetscFunctionBegin;
12659ac8732SStefano Zampini   if (mumps->schur_factored) {
12759ac8732SStefano Zampini     PetscFunctionReturn(0);
12859ac8732SStefano Zampini   }
12959ac8732SStefano Zampini   ierr = PetscBLASIntCast(mumps->id.size_schur,&B_N);CHKERRQ(ierr);
13059ac8732SStefano Zampini   ierr = PetscBLASIntCast(mumps->id.schur_lld,&B_slda);CHKERRQ(ierr);
13159ac8732SStefano Zampini   if (!mumps->sym) { /* MUMPS always return a full Schur matrix */
13259ac8732SStefano Zampini     if (!mumps->schur_pivots) {
13359ac8732SStefano Zampini       ierr = PetscMalloc1(B_N,&mumps->schur_pivots);CHKERRQ(ierr);
13459ac8732SStefano Zampini     }
13559ac8732SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
13659ac8732SStefano Zampini     PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,(PetscScalar*)mumps->id.schur,&B_slda,mumps->schur_pivots,&B_ierr));
13759ac8732SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
13859ac8732SStefano Zampini     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
13959ac8732SStefano Zampini   } else { /* either full or lower-triangular (not packed) */
14059ac8732SStefano Zampini     char ord[2];
14159ac8732SStefano Zampini     if (mumps->id.ICNTL(19) == 2 || mumps->id.ICNTL(19) == 3) { /* lower triangular stored by columns or full matrix */
14259ac8732SStefano Zampini       sprintf(ord,"L");
14359ac8732SStefano Zampini     } else { /* ICNTL(19) == 1 lower triangular stored by rows */
14459ac8732SStefano Zampini       sprintf(ord,"U");
14559ac8732SStefano Zampini     }
146a0b0af32SStefano Zampini     if (mumps->schur_sym == 2) {
14759ac8732SStefano Zampini       if (!mumps->schur_pivots) {
14859ac8732SStefano Zampini         PetscScalar  lwork;
14959ac8732SStefano Zampini 
15059ac8732SStefano Zampini         ierr = PetscMalloc1(B_N,&mumps->schur_pivots);CHKERRQ(ierr);
15159ac8732SStefano Zampini         mumps->schur_B_lwork=-1;
15259ac8732SStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
15359ac8732SStefano Zampini         PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_(ord,&B_N,(PetscScalar*)mumps->id.schur,&B_slda,mumps->schur_pivots,&lwork,&mumps->schur_B_lwork,&B_ierr));
15459ac8732SStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
15559ac8732SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYTRF Lapack routine %d",(int)B_ierr);
15659ac8732SStefano Zampini         ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&mumps->schur_B_lwork);CHKERRQ(ierr);
15759ac8732SStefano Zampini         ierr = PetscMalloc1(mumps->schur_B_lwork,&mumps->schur_work);CHKERRQ(ierr);
15859ac8732SStefano Zampini       }
15959ac8732SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
16059ac8732SStefano Zampini       PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_(ord,&B_N,(PetscScalar*)mumps->id.schur,&B_slda,mumps->schur_pivots,mumps->schur_work,&mumps->schur_B_lwork,&B_ierr));
16159ac8732SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
16259ac8732SStefano Zampini       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
16359ac8732SStefano Zampini     } else {
16459ac8732SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
16559ac8732SStefano Zampini       PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_(ord,&B_N,(PetscScalar*)mumps->id.schur,&B_slda,&B_ierr));
16659ac8732SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
16759ac8732SStefano Zampini       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
16859ac8732SStefano Zampini     }
16959ac8732SStefano Zampini   }
17059ac8732SStefano Zampini   mumps->schur_factored = PETSC_TRUE;
17159ac8732SStefano Zampini   PetscFunctionReturn(0);
17259ac8732SStefano Zampini }
17359ac8732SStefano Zampini 
17459ac8732SStefano Zampini static PetscErrorCode MatMumpsInvertSchur_Private(Mat_MUMPS* mumps)
17559ac8732SStefano Zampini {
17659ac8732SStefano Zampini   PetscBLASInt   B_N,B_ierr,B_slda;
17759ac8732SStefano Zampini   PetscErrorCode ierr;
17859ac8732SStefano Zampini 
17959ac8732SStefano Zampini   PetscFunctionBegin;
18059ac8732SStefano Zampini   ierr = MatMumpsFactorSchur_Private(mumps);CHKERRQ(ierr);
18159ac8732SStefano Zampini   ierr = PetscBLASIntCast(mumps->id.size_schur,&B_N);CHKERRQ(ierr);
18259ac8732SStefano Zampini   ierr = PetscBLASIntCast(mumps->id.schur_lld,&B_slda);CHKERRQ(ierr);
18359ac8732SStefano Zampini   if (!mumps->sym) { /* MUMPS always return a full Schur matrix */
18459ac8732SStefano Zampini     if (!mumps->schur_work) {
18559ac8732SStefano Zampini       PetscScalar lwork;
18659ac8732SStefano Zampini 
18759ac8732SStefano Zampini       mumps->schur_B_lwork = -1;
18859ac8732SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
18959ac8732SStefano Zampini       PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,(PetscScalar*)mumps->id.schur,&B_N,mumps->schur_pivots,&lwork,&mumps->schur_B_lwork,&B_ierr));
19059ac8732SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
19159ac8732SStefano Zampini       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
19259ac8732SStefano Zampini       ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&mumps->schur_B_lwork);CHKERRQ(ierr);
19359ac8732SStefano Zampini       ierr = PetscMalloc1(mumps->schur_B_lwork,&mumps->schur_work);CHKERRQ(ierr);
19459ac8732SStefano Zampini     }
19559ac8732SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
19659ac8732SStefano Zampini     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,(PetscScalar*)mumps->id.schur,&B_N,mumps->schur_pivots,mumps->schur_work,&mumps->schur_B_lwork,&B_ierr));
19759ac8732SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
19859ac8732SStefano Zampini     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
19959ac8732SStefano Zampini   } else { /* either full or lower-triangular (not packed) */
20059ac8732SStefano Zampini     char ord[2];
20159ac8732SStefano Zampini     if (mumps->id.ICNTL(19) == 2 || mumps->id.ICNTL(19) == 3) { /* lower triangular stored by columns or full matrix */
20259ac8732SStefano Zampini       sprintf(ord,"L");
20359ac8732SStefano Zampini     } else { /* ICNTL(19) == 1 lower triangular stored by rows */
20459ac8732SStefano Zampini       sprintf(ord,"U");
20559ac8732SStefano Zampini     }
206a0b0af32SStefano Zampini     if (mumps->schur_sym == 2) {
20759ac8732SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
20859ac8732SStefano Zampini       PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_(ord,&B_N,(PetscScalar*)mumps->id.schur,&B_N,mumps->schur_pivots,mumps->schur_work,&B_ierr));
20959ac8732SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
21059ac8732SStefano Zampini       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
21159ac8732SStefano Zampini     } else {
21259ac8732SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
21359ac8732SStefano Zampini       PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_(ord,&B_N,(PetscScalar*)mumps->id.schur,&B_N,&B_ierr));
21459ac8732SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
21559ac8732SStefano Zampini       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
21659ac8732SStefano Zampini     }
21759ac8732SStefano Zampini   }
21859ac8732SStefano Zampini   mumps->schur_inverted = PETSC_TRUE;
21959ac8732SStefano Zampini   PetscFunctionReturn(0);
22059ac8732SStefano Zampini }
22159ac8732SStefano Zampini 
222e807eca7SStefano Zampini static PetscErrorCode MatMumpsSolveSchur_Private(Mat_MUMPS* mumps, PetscBool sol_in_redrhs)
22359ac8732SStefano Zampini {
22459ac8732SStefano Zampini   PetscBLASInt   B_N,B_Nrhs,B_ierr,B_slda,B_rlda;
22559ac8732SStefano Zampini   PetscScalar    one=1.,zero=0.;
22659ac8732SStefano Zampini   PetscErrorCode ierr;
22759ac8732SStefano Zampini 
22859ac8732SStefano Zampini   PetscFunctionBegin;
22959ac8732SStefano Zampini   ierr = MatMumpsFactorSchur_Private(mumps);CHKERRQ(ierr);
230b5fa320bSStefano Zampini   ierr = PetscBLASIntCast(mumps->id.size_schur,&B_N);CHKERRQ(ierr);
231b5fa320bSStefano Zampini   ierr = PetscBLASIntCast(mumps->id.schur_lld,&B_slda);CHKERRQ(ierr);
232b5fa320bSStefano Zampini   ierr = PetscBLASIntCast(mumps->id.nrhs,&B_Nrhs);CHKERRQ(ierr);
233b5fa320bSStefano Zampini   ierr = PetscBLASIntCast(mumps->id.lredrhs,&B_rlda);CHKERRQ(ierr);
23459ac8732SStefano Zampini   if (mumps->schur_inverted) {
23559ac8732SStefano Zampini     PetscInt sizesol = B_Nrhs*B_N;
23659ac8732SStefano Zampini     if (!mumps->schur_sol || sizesol > mumps->schur_sizesol) {
23759ac8732SStefano Zampini       ierr = PetscFree(mumps->schur_sol);CHKERRQ(ierr);
23859ac8732SStefano Zampini       ierr = PetscMalloc1(sizesol,&mumps->schur_sol);CHKERRQ(ierr);
23959ac8732SStefano Zampini       mumps->schur_sizesol = sizesol;
240b5fa320bSStefano Zampini     }
24159ac8732SStefano Zampini     if (!mumps->sym) {
24259ac8732SStefano Zampini       char type[2];
243b5fa320bSStefano Zampini       if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
24459ac8732SStefano Zampini         if (!mumps->id.ICNTL(9)) { /* transpose solve */
24559ac8732SStefano Zampini           sprintf(type,"N");
246b5fa320bSStefano Zampini         } else {
24759ac8732SStefano Zampini           sprintf(type,"T");
248b5fa320bSStefano Zampini         }
24959ac8732SStefano Zampini       } else { /* stored by columns */
25059ac8732SStefano Zampini         if (!mumps->id.ICNTL(9)) { /* transpose solve */
25159ac8732SStefano Zampini           sprintf(type,"T");
25259ac8732SStefano Zampini         } else {
25359ac8732SStefano Zampini           sprintf(type,"N");
25459ac8732SStefano Zampini         }
25559ac8732SStefano Zampini       }
25659ac8732SStefano Zampini       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));
25759ac8732SStefano Zampini     } else {
25859ac8732SStefano Zampini       char ord[2];
25959ac8732SStefano Zampini       if (mumps->id.ICNTL(19) == 2 || mumps->id.ICNTL(19) == 3) { /* lower triangular stored by columns or full matrix */
26059ac8732SStefano Zampini         sprintf(ord,"L");
26159ac8732SStefano Zampini       } else { /* ICNTL(19) == 1 lower triangular stored by rows */
26259ac8732SStefano Zampini         sprintf(ord,"U");
26359ac8732SStefano Zampini       }
26459ac8732SStefano Zampini       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));
26559ac8732SStefano Zampini     }
266e807eca7SStefano Zampini     if (sol_in_redrhs) {
26759ac8732SStefano Zampini       ierr = PetscMemcpy(mumps->id.redrhs,mumps->schur_sol,sizesol*sizeof(PetscScalar));CHKERRQ(ierr);
268e807eca7SStefano Zampini     }
2693280fbe9SStefano Zampini   } else { /* Schur complement has not yet been inverted */
270a12f35bfSStefano Zampini     MumpsScalar *orhs=NULL;
271a12f35bfSStefano Zampini 
272a12f35bfSStefano Zampini     if (!sol_in_redrhs) {
273a12f35bfSStefano Zampini       PetscInt sizesol = B_Nrhs*B_N;
274a12f35bfSStefano Zampini       if (!mumps->schur_sol || sizesol > mumps->schur_sizesol) {
275a12f35bfSStefano Zampini         ierr = PetscFree(mumps->schur_sol);CHKERRQ(ierr);
276a12f35bfSStefano Zampini         ierr = PetscMalloc1(sizesol,&mumps->schur_sol);CHKERRQ(ierr);
277a12f35bfSStefano Zampini         mumps->schur_sizesol = sizesol;
278a12f35bfSStefano Zampini       }
279a12f35bfSStefano Zampini       orhs = mumps->id.redrhs;
280a12f35bfSStefano Zampini       ierr = PetscMemcpy(mumps->schur_sol,mumps->id.redrhs,sizesol*sizeof(PetscScalar));CHKERRQ(ierr);
281a12f35bfSStefano Zampini       mumps->id.redrhs = (MumpsScalar*)mumps->schur_sol;
282a12f35bfSStefano Zampini     }
28359ac8732SStefano Zampini     if (!mumps->sym) { /* MUMPS always return a full Schur matrix */
28459ac8732SStefano Zampini       char type[2];
28559ac8732SStefano Zampini       if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
28659ac8732SStefano Zampini         if (!mumps->id.ICNTL(9)) { /* transpose solve */
28759ac8732SStefano Zampini           sprintf(type,"N");
28859ac8732SStefano Zampini         } else {
28959ac8732SStefano Zampini           sprintf(type,"T");
29059ac8732SStefano Zampini         }
29159ac8732SStefano Zampini       } else { /* stored by columns */
29259ac8732SStefano Zampini         if (!mumps->id.ICNTL(9)) { /* transpose solve */
29359ac8732SStefano Zampini           sprintf(type,"T");
29459ac8732SStefano Zampini         } else {
29559ac8732SStefano Zampini           sprintf(type,"N");
29659ac8732SStefano Zampini         }
29759ac8732SStefano Zampini       }
29859ac8732SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
29959ac8732SStefano Zampini       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));
300b5fa320bSStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
301b5fa320bSStefano Zampini       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRS Lapack routine %d",(int)B_ierr);
302b5fa320bSStefano Zampini     } else { /* either full or lower-triangular (not packed) */
303b5fa320bSStefano Zampini       char ord[2];
304b5fa320bSStefano Zampini       if (mumps->id.ICNTL(19) == 2 || mumps->id.ICNTL(19) == 3) { /* lower triangular stored by columns or full matrix */
305b5fa320bSStefano Zampini         sprintf(ord,"L");
306b5fa320bSStefano Zampini       } else { /* ICNTL(19) == 1 lower triangular stored by rows */
307b5fa320bSStefano Zampini         sprintf(ord,"U");
308b5fa320bSStefano Zampini       }
309a0b0af32SStefano Zampini       if (mumps->schur_sym == 2) {
310b5fa320bSStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
31159ac8732SStefano Zampini         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));
312b5fa320bSStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
313b5fa320bSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRS Lapack routine %d",(int)B_ierr);
314b5fa320bSStefano Zampini       } else {
315b5fa320bSStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
31659ac8732SStefano Zampini         PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_(ord,&B_N,&B_Nrhs,(PetscScalar*)mumps->id.schur,&B_slda,(PetscScalar*)mumps->id.redrhs,&B_rlda,&B_ierr));
317b5fa320bSStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
318b5fa320bSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRS Lapack routine %d",(int)B_ierr);
319b5fa320bSStefano Zampini       }
320b5fa320bSStefano Zampini     }
321e807eca7SStefano Zampini     if (!sol_in_redrhs) {
322a12f35bfSStefano Zampini       mumps->id.redrhs = orhs;
323e807eca7SStefano Zampini     }
32459ac8732SStefano Zampini   }
325b5fa320bSStefano Zampini   PetscFunctionReturn(0);
326b5fa320bSStefano Zampini }
327b5fa320bSStefano Zampini 
328b8f61ee1SStefano Zampini static PetscErrorCode MatMumpsHandleSchur_Private(Mat_MUMPS* mumps, PetscBool expansion)
329b5fa320bSStefano Zampini {
330b5fa320bSStefano Zampini   PetscErrorCode ierr;
331b5fa320bSStefano Zampini 
332b5fa320bSStefano Zampini   PetscFunctionBegin;
333b5fa320bSStefano Zampini   if (!mumps->id.ICNTL(19)) { /* do nothing when Schur complement has not been computed */
334b5fa320bSStefano Zampini     PetscFunctionReturn(0);
335b5fa320bSStefano Zampini   }
336b8f61ee1SStefano Zampini   if (!expansion) { /* prepare for the condensation step */
337b5fa320bSStefano Zampini     PetscInt sizeredrhs = mumps->id.nrhs*mumps->id.size_schur;
338b5fa320bSStefano Zampini     /* allocate MUMPS internal array to store reduced right-hand sides */
339b5fa320bSStefano Zampini     if (!mumps->id.redrhs || sizeredrhs > mumps->sizeredrhs) {
340b5fa320bSStefano Zampini       ierr = PetscFree(mumps->id.redrhs);CHKERRQ(ierr);
341b5fa320bSStefano Zampini       mumps->id.lredrhs = mumps->id.size_schur;
342b5fa320bSStefano Zampini       ierr = PetscMalloc1(mumps->id.nrhs*mumps->id.lredrhs,&mumps->id.redrhs);CHKERRQ(ierr);
343b5fa320bSStefano Zampini       mumps->sizeredrhs = mumps->id.nrhs*mumps->id.lredrhs;
344b5fa320bSStefano Zampini     }
345b5fa320bSStefano Zampini     mumps->id.ICNTL(26) = 1; /* condensation phase */
346b5fa320bSStefano Zampini   } else { /* prepare for the expansion step */
347b8f61ee1SStefano Zampini     /* solve Schur complement (this has to be done by the MUMPS user, so basically us) */
348e807eca7SStefano Zampini     ierr = MatMumpsSolveSchur_Private(mumps,PETSC_TRUE);CHKERRQ(ierr);
349b5fa320bSStefano Zampini     mumps->id.ICNTL(26) = 2; /* expansion phase */
350b5fa320bSStefano Zampini     PetscMUMPS_c(&mumps->id);
351b5fa320bSStefano Zampini     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));
352b5fa320bSStefano Zampini     /* restore defaults */
353b5fa320bSStefano Zampini     mumps->id.ICNTL(26) = -1;
354d3d598ffSStefano Zampini     /* free MUMPS internal array for redrhs if we have solved for multiple rhs in order to save memory space */
355d3d598ffSStefano Zampini     if (mumps->id.nrhs > 1) {
356d3d598ffSStefano Zampini       ierr = PetscFree(mumps->id.redrhs);CHKERRQ(ierr);
357d3d598ffSStefano Zampini       mumps->id.lredrhs = 0;
358d3d598ffSStefano Zampini       mumps->sizeredrhs = 0;
359d3d598ffSStefano Zampini     }
360b5fa320bSStefano Zampini   }
361b5fa320bSStefano Zampini   PetscFunctionReturn(0);
362b5fa320bSStefano Zampini }
363b5fa320bSStefano Zampini 
364397b6df1SKris Buschelman /*
365d341cd04SHong Zhang   MatConvertToTriples_A_B - convert Petsc matrix to triples: row[nz], col[nz], val[nz]
366d341cd04SHong Zhang 
367397b6df1SKris Buschelman   input:
36867877ebaSShri Abhyankar     A       - matrix in aij,baij or sbaij (bs=1) format
369397b6df1SKris Buschelman     shift   - 0: C style output triple; 1: Fortran style output triple.
370bccb9932SShri Abhyankar     reuse   - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple
371bccb9932SShri Abhyankar               MAT_REUSE_MATRIX:   only the values in v array are updated
372397b6df1SKris Buschelman   output:
373397b6df1SKris Buschelman     nnz     - dim of r, c, and v (number of local nonzero entries of A)
374397b6df1SKris Buschelman     r, c, v - row and col index, matrix values (matrix triples)
375eb9baa12SBarry Smith 
376eb9baa12SBarry Smith   The returned values r, c, and sometimes v are obtained in a single PetscMalloc(). Then in MatDestroy_MUMPS() it is
377eb9baa12SBarry Smith   freed with PetscFree((mumps->irn);  This is not ideal code, the fact that v is ONLY sometimes part of mumps->irn means
378eb9baa12SBarry Smith   that the PetscMalloc() cannot easily be replaced with a PetscMalloc3().
379eb9baa12SBarry Smith 
380397b6df1SKris Buschelman  */
38116ebf90aSShri Abhyankar 
382bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
383b24902e0SBarry Smith {
384185f6596SHong Zhang   const PetscInt *ai,*aj,*ajj,M=A->rmap->n;
38567877ebaSShri Abhyankar   PetscInt       nz,rnz,i,j;
386dfbe8321SBarry Smith   PetscErrorCode ierr;
387c1490034SHong Zhang   PetscInt       *row,*col;
38816ebf90aSShri Abhyankar   Mat_SeqAIJ     *aa=(Mat_SeqAIJ*)A->data;
389397b6df1SKris Buschelman 
390397b6df1SKris Buschelman   PetscFunctionBegin;
39116ebf90aSShri Abhyankar   *v=aa->a;
392bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
3932205254eSKarl Rupp     nz   = aa->nz;
3942205254eSKarl Rupp     ai   = aa->i;
3952205254eSKarl Rupp     aj   = aa->j;
39616ebf90aSShri Abhyankar     *nnz = nz;
397785e854fSJed Brown     ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr);
398185f6596SHong Zhang     col  = row + nz;
399185f6596SHong Zhang 
40016ebf90aSShri Abhyankar     nz = 0;
40116ebf90aSShri Abhyankar     for (i=0; i<M; i++) {
40216ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
40367877ebaSShri Abhyankar       ajj = aj + ai[i];
40467877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
40567877ebaSShri Abhyankar         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
40616ebf90aSShri Abhyankar       }
40716ebf90aSShri Abhyankar     }
40816ebf90aSShri Abhyankar     *r = row; *c = col;
40916ebf90aSShri Abhyankar   }
41016ebf90aSShri Abhyankar   PetscFunctionReturn(0);
41116ebf90aSShri Abhyankar }
412397b6df1SKris Buschelman 
413bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
41467877ebaSShri Abhyankar {
41567877ebaSShri Abhyankar   Mat_SeqBAIJ    *aa=(Mat_SeqBAIJ*)A->data;
41633d57670SJed Brown   const PetscInt *ai,*aj,*ajj,bs2 = aa->bs2;
41733d57670SJed Brown   PetscInt       bs,M,nz,idx=0,rnz,i,j,k,m;
41867877ebaSShri Abhyankar   PetscErrorCode ierr;
41967877ebaSShri Abhyankar   PetscInt       *row,*col;
42067877ebaSShri Abhyankar 
42167877ebaSShri Abhyankar   PetscFunctionBegin;
42233d57670SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
42333d57670SJed Brown   M = A->rmap->N/bs;
424cf3759fdSShri Abhyankar   *v = aa->a;
425bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
426cf3759fdSShri Abhyankar     ai   = aa->i; aj = aa->j;
42767877ebaSShri Abhyankar     nz   = bs2*aa->nz;
42867877ebaSShri Abhyankar     *nnz = nz;
429785e854fSJed Brown     ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr);
430185f6596SHong Zhang     col  = row + nz;
431185f6596SHong Zhang 
43267877ebaSShri Abhyankar     for (i=0; i<M; i++) {
43367877ebaSShri Abhyankar       ajj = aj + ai[i];
43467877ebaSShri Abhyankar       rnz = ai[i+1] - ai[i];
43567877ebaSShri Abhyankar       for (k=0; k<rnz; k++) {
43667877ebaSShri Abhyankar         for (j=0; j<bs; j++) {
43767877ebaSShri Abhyankar           for (m=0; m<bs; m++) {
43867877ebaSShri Abhyankar             row[idx]   = i*bs + m + shift;
439cf3759fdSShri Abhyankar             col[idx++] = bs*(ajj[k]) + j + shift;
44067877ebaSShri Abhyankar           }
44167877ebaSShri Abhyankar         }
44267877ebaSShri Abhyankar       }
44367877ebaSShri Abhyankar     }
444cf3759fdSShri Abhyankar     *r = row; *c = col;
44567877ebaSShri Abhyankar   }
44667877ebaSShri Abhyankar   PetscFunctionReturn(0);
44767877ebaSShri Abhyankar }
44867877ebaSShri Abhyankar 
449bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
45016ebf90aSShri Abhyankar {
45167877ebaSShri Abhyankar   const PetscInt *ai, *aj,*ajj,M=A->rmap->n;
45267877ebaSShri Abhyankar   PetscInt       nz,rnz,i,j;
45316ebf90aSShri Abhyankar   PetscErrorCode ierr;
45416ebf90aSShri Abhyankar   PetscInt       *row,*col;
45516ebf90aSShri Abhyankar   Mat_SeqSBAIJ   *aa=(Mat_SeqSBAIJ*)A->data;
45616ebf90aSShri Abhyankar 
45716ebf90aSShri Abhyankar   PetscFunctionBegin;
458882afa5aSHong Zhang   *v = aa->a;
459bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
4602205254eSKarl Rupp     nz   = aa->nz;
4612205254eSKarl Rupp     ai   = aa->i;
4622205254eSKarl Rupp     aj   = aa->j;
4632205254eSKarl Rupp     *v   = aa->a;
46416ebf90aSShri Abhyankar     *nnz = nz;
465785e854fSJed Brown     ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr);
466185f6596SHong Zhang     col  = row + nz;
467185f6596SHong Zhang 
46816ebf90aSShri Abhyankar     nz = 0;
46916ebf90aSShri Abhyankar     for (i=0; i<M; i++) {
47016ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
47167877ebaSShri Abhyankar       ajj = aj + ai[i];
47267877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
47367877ebaSShri Abhyankar         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
47416ebf90aSShri Abhyankar       }
47516ebf90aSShri Abhyankar     }
47616ebf90aSShri Abhyankar     *r = row; *c = col;
47716ebf90aSShri Abhyankar   }
47816ebf90aSShri Abhyankar   PetscFunctionReturn(0);
47916ebf90aSShri Abhyankar }
48016ebf90aSShri Abhyankar 
481bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
48216ebf90aSShri Abhyankar {
48367877ebaSShri Abhyankar   const PetscInt    *ai,*aj,*ajj,*adiag,M=A->rmap->n;
48467877ebaSShri Abhyankar   PetscInt          nz,rnz,i,j;
48567877ebaSShri Abhyankar   const PetscScalar *av,*v1;
48616ebf90aSShri Abhyankar   PetscScalar       *val;
48716ebf90aSShri Abhyankar   PetscErrorCode    ierr;
48816ebf90aSShri Abhyankar   PetscInt          *row,*col;
489829b1710SHong Zhang   Mat_SeqAIJ        *aa=(Mat_SeqAIJ*)A->data;
49029b521d4Sstefano_zampini   PetscBool         missing;
49116ebf90aSShri Abhyankar 
49216ebf90aSShri Abhyankar   PetscFunctionBegin;
49316ebf90aSShri Abhyankar   ai   =aa->i; aj=aa->j;av=aa->a;
49416ebf90aSShri Abhyankar   adiag=aa->diag;
49529b521d4Sstefano_zampini   ierr  = MatMissingDiagonal_SeqAIJ(A,&missing,&i);CHKERRQ(ierr);
496bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
497829b1710SHong Zhang     /* count nz in the uppper triangular part of A */
498829b1710SHong Zhang     nz = 0;
49929b521d4Sstefano_zampini     if (missing) {
50029b521d4Sstefano_zampini       for (i=0; i<M; i++) {
50129b521d4Sstefano_zampini         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
50229b521d4Sstefano_zampini           for (j=ai[i];j<ai[i+1];j++) {
50329b521d4Sstefano_zampini             if (aj[j] < i) continue;
50429b521d4Sstefano_zampini             nz++;
50529b521d4Sstefano_zampini           }
50629b521d4Sstefano_zampini         } else {
50729b521d4Sstefano_zampini           nz += ai[i+1] - adiag[i];
50829b521d4Sstefano_zampini         }
50929b521d4Sstefano_zampini       }
51029b521d4Sstefano_zampini     } else {
511829b1710SHong Zhang       for (i=0; i<M; i++) nz += ai[i+1] - adiag[i];
51229b521d4Sstefano_zampini     }
51316ebf90aSShri Abhyankar     *nnz = nz;
514829b1710SHong Zhang 
515185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
516185f6596SHong Zhang     col  = row + nz;
517185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
518185f6596SHong Zhang 
51916ebf90aSShri Abhyankar     nz = 0;
52029b521d4Sstefano_zampini     if (missing) {
52129b521d4Sstefano_zampini       for (i=0; i<M; i++) {
52229b521d4Sstefano_zampini         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
52329b521d4Sstefano_zampini           for (j=ai[i];j<ai[i+1];j++) {
52429b521d4Sstefano_zampini             if (aj[j] < i) continue;
52529b521d4Sstefano_zampini             row[nz] = i+shift;
52629b521d4Sstefano_zampini             col[nz] = aj[j]+shift;
52729b521d4Sstefano_zampini             val[nz] = av[j];
52829b521d4Sstefano_zampini             nz++;
52929b521d4Sstefano_zampini           }
53029b521d4Sstefano_zampini         } else {
53129b521d4Sstefano_zampini           rnz = ai[i+1] - adiag[i];
53229b521d4Sstefano_zampini           ajj = aj + adiag[i];
53329b521d4Sstefano_zampini           v1  = av + adiag[i];
53429b521d4Sstefano_zampini           for (j=0; j<rnz; j++) {
53529b521d4Sstefano_zampini             row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
53629b521d4Sstefano_zampini           }
53729b521d4Sstefano_zampini         }
53829b521d4Sstefano_zampini       }
53929b521d4Sstefano_zampini     } else {
54016ebf90aSShri Abhyankar       for (i=0; i<M; i++) {
54116ebf90aSShri Abhyankar         rnz = ai[i+1] - adiag[i];
54267877ebaSShri Abhyankar         ajj = aj + adiag[i];
543cf3759fdSShri Abhyankar         v1  = av + adiag[i];
54467877ebaSShri Abhyankar         for (j=0; j<rnz; j++) {
54567877ebaSShri Abhyankar           row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
54616ebf90aSShri Abhyankar         }
54716ebf90aSShri Abhyankar       }
54829b521d4Sstefano_zampini     }
54916ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
550397b6df1SKris Buschelman   } else {
55116ebf90aSShri Abhyankar     nz = 0; val = *v;
55229b521d4Sstefano_zampini     if (missing) {
55316ebf90aSShri Abhyankar       for (i=0; i <M; i++) {
55429b521d4Sstefano_zampini         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
55529b521d4Sstefano_zampini           for (j=ai[i];j<ai[i+1];j++) {
55629b521d4Sstefano_zampini             if (aj[j] < i) continue;
55729b521d4Sstefano_zampini             val[nz++] = av[j];
55829b521d4Sstefano_zampini           }
55929b521d4Sstefano_zampini         } else {
56016ebf90aSShri Abhyankar           rnz = ai[i+1] - adiag[i];
56167877ebaSShri Abhyankar           v1  = av + adiag[i];
56267877ebaSShri Abhyankar           for (j=0; j<rnz; j++) {
56367877ebaSShri Abhyankar             val[nz++] = v1[j];
56416ebf90aSShri Abhyankar           }
56516ebf90aSShri Abhyankar         }
56616ebf90aSShri Abhyankar       }
56729b521d4Sstefano_zampini     } else {
56816ebf90aSShri Abhyankar       for (i=0; i <M; i++) {
56916ebf90aSShri Abhyankar         rnz = ai[i+1] - adiag[i];
57016ebf90aSShri Abhyankar         v1  = av + adiag[i];
57116ebf90aSShri Abhyankar         for (j=0; j<rnz; j++) {
57216ebf90aSShri Abhyankar           val[nz++] = v1[j];
57316ebf90aSShri Abhyankar         }
57416ebf90aSShri Abhyankar       }
57516ebf90aSShri Abhyankar     }
57629b521d4Sstefano_zampini   }
57716ebf90aSShri Abhyankar   PetscFunctionReturn(0);
57816ebf90aSShri Abhyankar }
57916ebf90aSShri Abhyankar 
580bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
58116ebf90aSShri Abhyankar {
58216ebf90aSShri Abhyankar   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
58316ebf90aSShri Abhyankar   PetscErrorCode    ierr;
58416ebf90aSShri Abhyankar   PetscInt          rstart,nz,i,j,jj,irow,countA,countB;
58516ebf90aSShri Abhyankar   PetscInt          *row,*col;
58616ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
58716ebf90aSShri Abhyankar   PetscScalar       *val;
588397b6df1SKris Buschelman   Mat_MPISBAIJ      *mat = (Mat_MPISBAIJ*)A->data;
589397b6df1SKris Buschelman   Mat_SeqSBAIJ      *aa  = (Mat_SeqSBAIJ*)(mat->A)->data;
590397b6df1SKris Buschelman   Mat_SeqBAIJ       *bb  = (Mat_SeqBAIJ*)(mat->B)->data;
59116ebf90aSShri Abhyankar 
59216ebf90aSShri Abhyankar   PetscFunctionBegin;
593d0f46423SBarry Smith   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
594397b6df1SKris Buschelman   av=aa->a; bv=bb->a;
595397b6df1SKris Buschelman 
5962205254eSKarl Rupp   garray = mat->garray;
5972205254eSKarl Rupp 
598bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
59916ebf90aSShri Abhyankar     nz   = aa->nz + bb->nz;
60016ebf90aSShri Abhyankar     *nnz = nz;
601185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
602185f6596SHong Zhang     col  = row + nz;
603185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
604185f6596SHong Zhang 
605397b6df1SKris Buschelman     *r = row; *c = col; *v = val;
606397b6df1SKris Buschelman   } else {
607397b6df1SKris Buschelman     row = *r; col = *c; val = *v;
608397b6df1SKris Buschelman   }
609397b6df1SKris Buschelman 
610028e57e8SHong Zhang   jj = 0; irow = rstart;
611397b6df1SKris Buschelman   for (i=0; i<m; i++) {
612397b6df1SKris Buschelman     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
613397b6df1SKris Buschelman     countA = ai[i+1] - ai[i];
614397b6df1SKris Buschelman     countB = bi[i+1] - bi[i];
615397b6df1SKris Buschelman     bjj    = bj + bi[i];
61616ebf90aSShri Abhyankar     v1     = av + ai[i];
61716ebf90aSShri Abhyankar     v2     = bv + bi[i];
618397b6df1SKris Buschelman 
619397b6df1SKris Buschelman     /* A-part */
620397b6df1SKris Buschelman     for (j=0; j<countA; j++) {
621bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
622397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
623397b6df1SKris Buschelman       }
62416ebf90aSShri Abhyankar       val[jj++] = v1[j];
625397b6df1SKris Buschelman     }
62616ebf90aSShri Abhyankar 
62716ebf90aSShri Abhyankar     /* B-part */
62816ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
629bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
630397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
631397b6df1SKris Buschelman       }
63216ebf90aSShri Abhyankar       val[jj++] = v2[j];
63316ebf90aSShri Abhyankar     }
63416ebf90aSShri Abhyankar     irow++;
63516ebf90aSShri Abhyankar   }
63616ebf90aSShri Abhyankar   PetscFunctionReturn(0);
63716ebf90aSShri Abhyankar }
63816ebf90aSShri Abhyankar 
639bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
64016ebf90aSShri Abhyankar {
64116ebf90aSShri Abhyankar   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
64216ebf90aSShri Abhyankar   PetscErrorCode    ierr;
64316ebf90aSShri Abhyankar   PetscInt          rstart,nz,i,j,jj,irow,countA,countB;
64416ebf90aSShri Abhyankar   PetscInt          *row,*col;
64516ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
64616ebf90aSShri Abhyankar   PetscScalar       *val;
64716ebf90aSShri Abhyankar   Mat_MPIAIJ        *mat = (Mat_MPIAIJ*)A->data;
64816ebf90aSShri Abhyankar   Mat_SeqAIJ        *aa  = (Mat_SeqAIJ*)(mat->A)->data;
64916ebf90aSShri Abhyankar   Mat_SeqAIJ        *bb  = (Mat_SeqAIJ*)(mat->B)->data;
65016ebf90aSShri Abhyankar 
65116ebf90aSShri Abhyankar   PetscFunctionBegin;
65216ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
65316ebf90aSShri Abhyankar   av=aa->a; bv=bb->a;
65416ebf90aSShri Abhyankar 
6552205254eSKarl Rupp   garray = mat->garray;
6562205254eSKarl Rupp 
657bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
65816ebf90aSShri Abhyankar     nz   = aa->nz + bb->nz;
65916ebf90aSShri Abhyankar     *nnz = nz;
660185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
661185f6596SHong Zhang     col  = row + nz;
662185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
663185f6596SHong Zhang 
66416ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
66516ebf90aSShri Abhyankar   } else {
66616ebf90aSShri Abhyankar     row = *r; col = *c; val = *v;
66716ebf90aSShri Abhyankar   }
66816ebf90aSShri Abhyankar 
66916ebf90aSShri Abhyankar   jj = 0; irow = rstart;
67016ebf90aSShri Abhyankar   for (i=0; i<m; i++) {
67116ebf90aSShri Abhyankar     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
67216ebf90aSShri Abhyankar     countA = ai[i+1] - ai[i];
67316ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
67416ebf90aSShri Abhyankar     bjj    = bj + bi[i];
67516ebf90aSShri Abhyankar     v1     = av + ai[i];
67616ebf90aSShri Abhyankar     v2     = bv + bi[i];
67716ebf90aSShri Abhyankar 
67816ebf90aSShri Abhyankar     /* A-part */
67916ebf90aSShri Abhyankar     for (j=0; j<countA; j++) {
680bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
68116ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
68216ebf90aSShri Abhyankar       }
68316ebf90aSShri Abhyankar       val[jj++] = v1[j];
68416ebf90aSShri Abhyankar     }
68516ebf90aSShri Abhyankar 
68616ebf90aSShri Abhyankar     /* B-part */
68716ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
688bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
68916ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
69016ebf90aSShri Abhyankar       }
69116ebf90aSShri Abhyankar       val[jj++] = v2[j];
69216ebf90aSShri Abhyankar     }
69316ebf90aSShri Abhyankar     irow++;
69416ebf90aSShri Abhyankar   }
69516ebf90aSShri Abhyankar   PetscFunctionReturn(0);
69616ebf90aSShri Abhyankar }
69716ebf90aSShri Abhyankar 
698bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
69967877ebaSShri Abhyankar {
70067877ebaSShri Abhyankar   Mat_MPIBAIJ       *mat    = (Mat_MPIBAIJ*)A->data;
70167877ebaSShri Abhyankar   Mat_SeqBAIJ       *aa     = (Mat_SeqBAIJ*)(mat->A)->data;
70267877ebaSShri Abhyankar   Mat_SeqBAIJ       *bb     = (Mat_SeqBAIJ*)(mat->B)->data;
70367877ebaSShri Abhyankar   const PetscInt    *ai     = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj;
704d985c460SShri Abhyankar   const PetscInt    *garray = mat->garray,mbs=mat->mbs,rstart=A->rmap->rstart;
70533d57670SJed Brown   const PetscInt    bs2=mat->bs2;
70667877ebaSShri Abhyankar   PetscErrorCode    ierr;
70733d57670SJed Brown   PetscInt          bs,nz,i,j,k,n,jj,irow,countA,countB,idx;
70867877ebaSShri Abhyankar   PetscInt          *row,*col;
70967877ebaSShri Abhyankar   const PetscScalar *av=aa->a, *bv=bb->a,*v1,*v2;
71067877ebaSShri Abhyankar   PetscScalar       *val;
71167877ebaSShri Abhyankar 
71267877ebaSShri Abhyankar   PetscFunctionBegin;
71333d57670SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
714bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
71567877ebaSShri Abhyankar     nz   = bs2*(aa->nz + bb->nz);
71667877ebaSShri Abhyankar     *nnz = nz;
717185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
718185f6596SHong Zhang     col  = row + nz;
719185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
720185f6596SHong Zhang 
72167877ebaSShri Abhyankar     *r = row; *c = col; *v = val;
72267877ebaSShri Abhyankar   } else {
72367877ebaSShri Abhyankar     row = *r; col = *c; val = *v;
72467877ebaSShri Abhyankar   }
72567877ebaSShri Abhyankar 
726d985c460SShri Abhyankar   jj = 0; irow = rstart;
72767877ebaSShri Abhyankar   for (i=0; i<mbs; i++) {
72867877ebaSShri Abhyankar     countA = ai[i+1] - ai[i];
72967877ebaSShri Abhyankar     countB = bi[i+1] - bi[i];
73067877ebaSShri Abhyankar     ajj    = aj + ai[i];
73167877ebaSShri Abhyankar     bjj    = bj + bi[i];
73267877ebaSShri Abhyankar     v1     = av + bs2*ai[i];
73367877ebaSShri Abhyankar     v2     = bv + bs2*bi[i];
73467877ebaSShri Abhyankar 
73567877ebaSShri Abhyankar     idx = 0;
73667877ebaSShri Abhyankar     /* A-part */
73767877ebaSShri Abhyankar     for (k=0; k<countA; k++) {
73867877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
73967877ebaSShri Abhyankar         for (n=0; n<bs; n++) {
740bccb9932SShri Abhyankar           if (reuse == MAT_INITIAL_MATRIX) {
741d985c460SShri Abhyankar             row[jj] = irow + n + shift;
742d985c460SShri Abhyankar             col[jj] = rstart + bs*ajj[k] + j + shift;
74367877ebaSShri Abhyankar           }
74467877ebaSShri Abhyankar           val[jj++] = v1[idx++];
74567877ebaSShri Abhyankar         }
74667877ebaSShri Abhyankar       }
74767877ebaSShri Abhyankar     }
74867877ebaSShri Abhyankar 
74967877ebaSShri Abhyankar     idx = 0;
75067877ebaSShri Abhyankar     /* B-part */
75167877ebaSShri Abhyankar     for (k=0; k<countB; k++) {
75267877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
75367877ebaSShri Abhyankar         for (n=0; n<bs; n++) {
754bccb9932SShri Abhyankar           if (reuse == MAT_INITIAL_MATRIX) {
755d985c460SShri Abhyankar             row[jj] = irow + n + shift;
756d985c460SShri Abhyankar             col[jj] = bs*garray[bjj[k]] + j + shift;
75767877ebaSShri Abhyankar           }
758d985c460SShri Abhyankar           val[jj++] = v2[idx++];
75967877ebaSShri Abhyankar         }
76067877ebaSShri Abhyankar       }
76167877ebaSShri Abhyankar     }
762d985c460SShri Abhyankar     irow += bs;
76367877ebaSShri Abhyankar   }
76467877ebaSShri Abhyankar   PetscFunctionReturn(0);
76567877ebaSShri Abhyankar }
76667877ebaSShri Abhyankar 
767bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
76816ebf90aSShri Abhyankar {
76916ebf90aSShri Abhyankar   const PetscInt    *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
77016ebf90aSShri Abhyankar   PetscErrorCode    ierr;
771e0bace9bSHong Zhang   PetscInt          rstart,nz,nza,nzb,i,j,jj,irow,countA,countB;
77216ebf90aSShri Abhyankar   PetscInt          *row,*col;
77316ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
77416ebf90aSShri Abhyankar   PetscScalar       *val;
77516ebf90aSShri Abhyankar   Mat_MPIAIJ        *mat =  (Mat_MPIAIJ*)A->data;
77616ebf90aSShri Abhyankar   Mat_SeqAIJ        *aa  =(Mat_SeqAIJ*)(mat->A)->data;
77716ebf90aSShri Abhyankar   Mat_SeqAIJ        *bb  =(Mat_SeqAIJ*)(mat->B)->data;
77816ebf90aSShri Abhyankar 
77916ebf90aSShri Abhyankar   PetscFunctionBegin;
78016ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j; adiag=aa->diag;
78116ebf90aSShri Abhyankar   bi=bb->i; bj=bb->j; garray = mat->garray;
78216ebf90aSShri Abhyankar   av=aa->a; bv=bb->a;
7832205254eSKarl Rupp 
78416ebf90aSShri Abhyankar   rstart = A->rmap->rstart;
78516ebf90aSShri Abhyankar 
786bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
787e0bace9bSHong Zhang     nza = 0;    /* num of upper triangular entries in mat->A, including diagonals */
788e0bace9bSHong Zhang     nzb = 0;    /* num of upper triangular entries in mat->B */
78916ebf90aSShri Abhyankar     for (i=0; i<m; i++) {
790e0bace9bSHong Zhang       nza   += (ai[i+1] - adiag[i]);
79116ebf90aSShri Abhyankar       countB = bi[i+1] - bi[i];
79216ebf90aSShri Abhyankar       bjj    = bj + bi[i];
793e0bace9bSHong Zhang       for (j=0; j<countB; j++) {
794e0bace9bSHong Zhang         if (garray[bjj[j]] > rstart) nzb++;
795e0bace9bSHong Zhang       }
796e0bace9bSHong Zhang     }
79716ebf90aSShri Abhyankar 
798e0bace9bSHong Zhang     nz   = nza + nzb; /* total nz of upper triangular part of mat */
79916ebf90aSShri Abhyankar     *nnz = nz;
800185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
801185f6596SHong Zhang     col  = row + nz;
802185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
803185f6596SHong Zhang 
80416ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
80516ebf90aSShri Abhyankar   } else {
80616ebf90aSShri Abhyankar     row = *r; col = *c; val = *v;
80716ebf90aSShri Abhyankar   }
80816ebf90aSShri Abhyankar 
80916ebf90aSShri Abhyankar   jj = 0; irow = rstart;
81016ebf90aSShri Abhyankar   for (i=0; i<m; i++) {
81116ebf90aSShri Abhyankar     ajj    = aj + adiag[i];                 /* ptr to the beginning of the diagonal of this row */
81216ebf90aSShri Abhyankar     v1     = av + adiag[i];
81316ebf90aSShri Abhyankar     countA = ai[i+1] - adiag[i];
81416ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
81516ebf90aSShri Abhyankar     bjj    = bj + bi[i];
81616ebf90aSShri Abhyankar     v2     = bv + bi[i];
81716ebf90aSShri Abhyankar 
81816ebf90aSShri Abhyankar     /* A-part */
81916ebf90aSShri Abhyankar     for (j=0; j<countA; j++) {
820bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
82116ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
82216ebf90aSShri Abhyankar       }
82316ebf90aSShri Abhyankar       val[jj++] = v1[j];
82416ebf90aSShri Abhyankar     }
82516ebf90aSShri Abhyankar 
82616ebf90aSShri Abhyankar     /* B-part */
82716ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
82816ebf90aSShri Abhyankar       if (garray[bjj[j]] > rstart) {
829bccb9932SShri Abhyankar         if (reuse == MAT_INITIAL_MATRIX) {
83016ebf90aSShri Abhyankar           row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
83116ebf90aSShri Abhyankar         }
83216ebf90aSShri Abhyankar         val[jj++] = v2[j];
83316ebf90aSShri Abhyankar       }
834397b6df1SKris Buschelman     }
835397b6df1SKris Buschelman     irow++;
836397b6df1SKris Buschelman   }
837397b6df1SKris Buschelman   PetscFunctionReturn(0);
838397b6df1SKris Buschelman }
839397b6df1SKris Buschelman 
840dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A)
841dfbe8321SBarry Smith {
842e69c285eSBarry Smith   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->data;
843dfbe8321SBarry Smith   PetscErrorCode ierr;
844b24902e0SBarry Smith 
845397b6df1SKris Buschelman   PetscFunctionBegin;
846a5e57a09SHong Zhang   ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr);
847a5e57a09SHong Zhang   ierr = VecScatterDestroy(&mumps->scat_rhs);CHKERRQ(ierr);
848a5e57a09SHong Zhang   ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
849801fbe65SHong Zhang   ierr = VecDestroy(&mumps->b_seq);CHKERRQ(ierr);
850a5e57a09SHong Zhang   ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr);
851a5e57a09SHong Zhang   ierr = PetscFree(mumps->id.perm_in);CHKERRQ(ierr);
852a5e57a09SHong Zhang   ierr = PetscFree(mumps->irn);CHKERRQ(ierr);
853b34f08ffSHong Zhang   ierr = PetscFree(mumps->info);CHKERRQ(ierr);
85459ac8732SStefano Zampini   ierr = MatMumpsResetSchur_Private(mumps);CHKERRQ(ierr);
855a5e57a09SHong Zhang   mumps->id.job = JOB_END;
856a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
8576f3cc6f9SBarry Smith   ierr = MPI_Comm_free(&mumps->comm_mumps);CHKERRQ(ierr);
858e69c285eSBarry Smith   ierr = PetscFree(A->data);CHKERRQ(ierr);
859bf0cc555SLisandro Dalcin 
86097969023SHong Zhang   /* clear composed functions */
861bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr);
8625a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorSetSchurIS_C",NULL);CHKERRQ(ierr);
8635a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorInvertSchurComplement_C",NULL);CHKERRQ(ierr);
8645a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorCreateSchurComplement_C",NULL);CHKERRQ(ierr);
8655a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSchurComplement_C",NULL);CHKERRQ(ierr);
8665a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorSolveSchurComplement_C",NULL);CHKERRQ(ierr);
8675a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorSolveSchurComplementTranspose_C",NULL);CHKERRQ(ierr);
8686dba178dSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorFactorizeSchurComplement_C",NULL);CHKERRQ(ierr);
869e8ade678SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorSetSchurComplementSolverType_C",NULL);CHKERRQ(ierr);
870bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetIcntl_C",NULL);CHKERRQ(ierr);
871bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetIcntl_C",NULL);CHKERRQ(ierr);
872bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetCntl_C",NULL);CHKERRQ(ierr);
873bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetCntl_C",NULL);CHKERRQ(ierr);
874ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfo_C",NULL);CHKERRQ(ierr);
875ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfog_C",NULL);CHKERRQ(ierr);
876ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfo_C",NULL);CHKERRQ(ierr);
877ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfog_C",NULL);CHKERRQ(ierr);
878397b6df1SKris Buschelman   PetscFunctionReturn(0);
879397b6df1SKris Buschelman }
880397b6df1SKris Buschelman 
881b24902e0SBarry Smith PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x)
882b24902e0SBarry Smith {
883e69c285eSBarry Smith   Mat_MUMPS        *mumps=(Mat_MUMPS*)A->data;
884d54de34fSKris Buschelman   PetscScalar      *array;
88567877ebaSShri Abhyankar   Vec              b_seq;
886329ec9b3SHong Zhang   IS               is_iden,is_petsc;
887dfbe8321SBarry Smith   PetscErrorCode   ierr;
888329ec9b3SHong Zhang   PetscInt         i;
889cc86f929SStefano Zampini   PetscBool        second_solve = PETSC_FALSE;
890883f2eb9SBarry Smith   static PetscBool cite1 = PETSC_FALSE,cite2 = PETSC_FALSE;
891397b6df1SKris Buschelman 
892397b6df1SKris Buschelman   PetscFunctionBegin;
893883f2eb9SBarry Smith   ierr = 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);CHKERRQ(ierr);
894883f2eb9SBarry Smith   ierr = 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);CHKERRQ(ierr);
8952aca8efcSHong Zhang 
896603e8f96SBarry Smith   if (A->factorerrortype) {
8972aca8efcSHong Zhang     ierr = PetscInfo2(A,"MatSolve is called with singular matrix factor, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
8982aca8efcSHong Zhang     ierr = VecSetInf(x);CHKERRQ(ierr);
8992aca8efcSHong Zhang     PetscFunctionReturn(0);
9002aca8efcSHong Zhang   }
9012aca8efcSHong Zhang 
902a5e57a09SHong Zhang   mumps->id.nrhs = 1;
903a5e57a09SHong Zhang   b_seq          = mumps->b_seq;
904a5e57a09SHong Zhang   if (mumps->size > 1) {
905329ec9b3SHong Zhang     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
906a5e57a09SHong Zhang     ierr = VecScatterBegin(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
907a5e57a09SHong Zhang     ierr = VecScatterEnd(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
908a5e57a09SHong Zhang     if (!mumps->myid) {ierr = VecGetArray(b_seq,&array);CHKERRQ(ierr);}
909397b6df1SKris Buschelman   } else {  /* size == 1 */
910397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
911397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
912397b6df1SKris Buschelman   }
913a5e57a09SHong Zhang   if (!mumps->myid) { /* define rhs on the host */
914a5e57a09SHong Zhang     mumps->id.nrhs = 1;
915940cd9d6SSatish Balay     mumps->id.rhs = (MumpsScalar*)array;
916397b6df1SKris Buschelman   }
917397b6df1SKris Buschelman 
918cc86f929SStefano Zampini   /*
919cc86f929SStefano Zampini      handle condensation step of Schur complement (if any)
920cc86f929SStefano Zampini      We set by default ICNTL(26) == -1 when Schur indices have been provided by the user.
921cc86f929SStefano Zampini      According to MUMPS (5.0.0) manual, any value should be harmful during the factorization phase
922cc86f929SStefano Zampini      Unless the user provides a valid value for ICNTL(26), MatSolve and MatMatSolve routines solve the full system.
923cc86f929SStefano Zampini      This requires an extra call to PetscMUMPS_c and the computation of the factors for S
924cc86f929SStefano Zampini   */
925cc86f929SStefano Zampini   if (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2) {
926241dbb5eSStefano Zampini     if (mumps->size > 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Parallel Schur complements not yet supported from PETSc\n");
927cc86f929SStefano Zampini     second_solve = PETSC_TRUE;
928b8f61ee1SStefano Zampini     ierr = MatMumpsHandleSchur_Private(mumps,PETSC_FALSE);CHKERRQ(ierr);
929cc86f929SStefano Zampini   }
930397b6df1SKris Buschelman   /* solve phase */
931329ec9b3SHong Zhang   /*-------------*/
932a5e57a09SHong Zhang   mumps->id.job = JOB_SOLVE;
933a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
934a5e57a09SHong Zhang   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));
935397b6df1SKris Buschelman 
936b5fa320bSStefano Zampini   /* handle expansion step of Schur complement (if any) */
937cc86f929SStefano Zampini   if (second_solve) {
938b8f61ee1SStefano Zampini     ierr = MatMumpsHandleSchur_Private(mumps,PETSC_TRUE);CHKERRQ(ierr);
939cc86f929SStefano Zampini   }
940b5fa320bSStefano Zampini 
941a5e57a09SHong Zhang   if (mumps->size > 1) { /* convert mumps distributed solution to petsc mpi x */
942a5e57a09SHong Zhang     if (mumps->scat_sol && mumps->ICNTL9_pre != mumps->id.ICNTL(9)) {
943a5e57a09SHong Zhang       /* when id.ICNTL(9) changes, the contents of lsol_loc may change (not its size, lsol_loc), recreates scat_sol */
944a5e57a09SHong Zhang       ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
945397b6df1SKris Buschelman     }
946a5e57a09SHong Zhang     if (!mumps->scat_sol) { /* create scatter scat_sol */
947a5e57a09SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,mumps->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */
948a5e57a09SHong Zhang       for (i=0; i<mumps->id.lsol_loc; i++) {
949a5e57a09SHong Zhang         mumps->id.isol_loc[i] -= 1; /* change Fortran style to C style */
950a5e57a09SHong Zhang       }
951a5e57a09SHong Zhang       ierr = ISCreateGeneral(PETSC_COMM_SELF,mumps->id.lsol_loc,mumps->id.isol_loc,PETSC_COPY_VALUES,&is_petsc);CHKERRQ(ierr);  /* to */
952a5e57a09SHong Zhang       ierr = VecScatterCreate(mumps->x_seq,is_iden,x,is_petsc,&mumps->scat_sol);CHKERRQ(ierr);
9536bf464f9SBarry Smith       ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
9546bf464f9SBarry Smith       ierr = ISDestroy(&is_petsc);CHKERRQ(ierr);
9552205254eSKarl Rupp 
956a5e57a09SHong Zhang       mumps->ICNTL9_pre = mumps->id.ICNTL(9); /* save current value of id.ICNTL(9) */
957397b6df1SKris Buschelman     }
958a5e57a09SHong Zhang 
959a5e57a09SHong Zhang     ierr = VecScatterBegin(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
960a5e57a09SHong Zhang     ierr = VecScatterEnd(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
961329ec9b3SHong Zhang   }
962397b6df1SKris Buschelman   PetscFunctionReturn(0);
963397b6df1SKris Buschelman }
964397b6df1SKris Buschelman 
96551d5961aSHong Zhang PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x)
96651d5961aSHong Zhang {
967e69c285eSBarry Smith   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->data;
96851d5961aSHong Zhang   PetscErrorCode ierr;
96951d5961aSHong Zhang 
97051d5961aSHong Zhang   PetscFunctionBegin;
971a5e57a09SHong Zhang   mumps->id.ICNTL(9) = 0;
9720ad0caddSJed Brown   ierr = MatSolve_MUMPS(A,b,x);CHKERRQ(ierr);
973a5e57a09SHong Zhang   mumps->id.ICNTL(9) = 1;
97451d5961aSHong Zhang   PetscFunctionReturn(0);
97551d5961aSHong Zhang }
97651d5961aSHong Zhang 
977e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X)
978e0b74bf9SHong Zhang {
979bda8bf91SBarry Smith   PetscErrorCode ierr;
980bda8bf91SBarry Smith   PetscBool      flg;
981e69c285eSBarry Smith   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->data;
982334c5f61SHong Zhang   PetscInt       i,nrhs,M;
9832cd7d884SHong Zhang   PetscScalar    *array,*bray;
984bda8bf91SBarry Smith 
985e0b74bf9SHong Zhang   PetscFunctionBegin;
9860298fd71SBarry Smith   ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
987801fbe65SHong Zhang   if (!flg) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
9880298fd71SBarry Smith   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
989801fbe65SHong Zhang   if (!flg) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
990801fbe65SHong Zhang   if (B->rmap->n != X->rmap->n) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B and X must have same row distribution");
9914e34a73bSHong Zhang 
9922cd7d884SHong Zhang   ierr = MatGetSize(B,&M,&nrhs);CHKERRQ(ierr);
993334c5f61SHong Zhang   mumps->id.nrhs = nrhs;
994334c5f61SHong Zhang   mumps->id.lrhs = M;
9954e34a73bSHong Zhang 
9962cd7d884SHong Zhang   if (mumps->size == 1) {
997e94cce23SStefano Zampini     PetscBool second_solve = PETSC_FALSE;
9982cd7d884SHong Zhang     /* copy B to X */
9992cd7d884SHong Zhang     ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr);
10002cd7d884SHong Zhang     ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr);
10016444a565SStefano Zampini     ierr = PetscMemcpy(array,bray,M*nrhs*sizeof(PetscScalar));CHKERRQ(ierr);
10022cd7d884SHong Zhang     ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr);
1003940cd9d6SSatish Balay     mumps->id.rhs = (MumpsScalar*)array;
1004801fbe65SHong Zhang 
1005e94cce23SStefano Zampini     /* handle condensation step of Schur complement (if any) */
1006e94cce23SStefano Zampini     if (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2) {
1007e94cce23SStefano Zampini       second_solve = PETSC_TRUE;
1008e94cce23SStefano Zampini       ierr = MatMumpsHandleSchur_Private(mumps,PETSC_FALSE);CHKERRQ(ierr);
1009e94cce23SStefano Zampini     }
10102cd7d884SHong Zhang     /* solve phase */
10112cd7d884SHong Zhang     /*-------------*/
10122cd7d884SHong Zhang     mumps->id.job = JOB_SOLVE;
10132cd7d884SHong Zhang     PetscMUMPS_c(&mumps->id);
10142cd7d884SHong Zhang     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));
1015b5fa320bSStefano Zampini 
1016b5fa320bSStefano Zampini     /* handle expansion step of Schur complement (if any) */
1017e94cce23SStefano Zampini     if (second_solve) {
1018b8f61ee1SStefano Zampini       ierr = MatMumpsHandleSchur_Private(mumps,PETSC_TRUE);CHKERRQ(ierr);
1019e94cce23SStefano Zampini     }
10202cd7d884SHong Zhang     ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr);
1021334c5f61SHong Zhang   } else {  /*--------- parallel case --------*/
102271aed81dSHong Zhang     PetscInt       lsol_loc,nlsol_loc,*isol_loc,*idx,*iidx,*idxx,*isol_loc_save;
10231070efccSSatish Balay     MumpsScalar    *sol_loc,*sol_loc_save;
1024801fbe65SHong Zhang     IS             is_to,is_from;
1025334c5f61SHong Zhang     PetscInt       k,proc,j,m;
1026801fbe65SHong Zhang     const PetscInt *rstart;
1027334c5f61SHong Zhang     Vec            v_mpi,b_seq,x_seq;
1028334c5f61SHong Zhang     VecScatter     scat_rhs,scat_sol;
1029801fbe65SHong Zhang 
1030*38be02acSStefano Zampini     if (mumps->size > 1 && mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Parallel Schur complements not yet supported from PETSc\n");
1031241dbb5eSStefano Zampini 
1032801fbe65SHong Zhang     /* create x_seq to hold local solution */
103371aed81dSHong Zhang     isol_loc_save = mumps->id.isol_loc; /* save it for MatSovle() */
103471aed81dSHong Zhang     sol_loc_save  = mumps->id.sol_loc;
1035801fbe65SHong Zhang 
103671aed81dSHong Zhang     lsol_loc  = mumps->id.INFO(23);
103771aed81dSHong Zhang     nlsol_loc = nrhs*lsol_loc;     /* length of sol_loc */
103871aed81dSHong Zhang     ierr = PetscMalloc2(nlsol_loc,&sol_loc,nlsol_loc,&isol_loc);CHKERRQ(ierr);
1039940cd9d6SSatish Balay     mumps->id.sol_loc = (MumpsScalar*)sol_loc;
1040801fbe65SHong Zhang     mumps->id.isol_loc = isol_loc;
1041801fbe65SHong Zhang 
10421070efccSSatish Balay     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,nlsol_loc,(PetscScalar*)sol_loc,&x_seq);CHKERRQ(ierr);
10432cd7d884SHong Zhang 
104474f0fcc7SHong Zhang     /* copy rhs matrix B into vector v_mpi */
1045334c5f61SHong Zhang     ierr = MatGetLocalSize(B,&m,NULL);CHKERRQ(ierr);
1046801fbe65SHong Zhang     ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr);
104774f0fcc7SHong Zhang     ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)B),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi);CHKERRQ(ierr);
1048801fbe65SHong Zhang     ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr);
1049801fbe65SHong Zhang 
1050334c5f61SHong Zhang     /* scatter v_mpi to b_seq because MUMPS only supports centralized rhs */
105174f0fcc7SHong Zhang     /* idx: maps from k-th index of v_mpi to (i,j)-th global entry of B;
1052801fbe65SHong Zhang       iidx: inverse of idx, will be used by scattering xx_seq -> X       */
1053801fbe65SHong Zhang     ierr = PetscMalloc2(nrhs*M,&idx,nrhs*M,&iidx);CHKERRQ(ierr);
1054801fbe65SHong Zhang     ierr = MatGetOwnershipRanges(B,&rstart);CHKERRQ(ierr);
1055801fbe65SHong Zhang     k = 0;
1056801fbe65SHong Zhang     for (proc=0; proc<mumps->size; proc++){
1057801fbe65SHong Zhang       for (j=0; j<nrhs; j++){
1058801fbe65SHong Zhang         for (i=rstart[proc]; i<rstart[proc+1]; i++){
1059801fbe65SHong Zhang           iidx[j*M + i] = k;
1060801fbe65SHong Zhang           idx[k++]      = j*M + i;
1061801fbe65SHong Zhang         }
1062801fbe65SHong Zhang       }
10632cd7d884SHong Zhang     }
10642cd7d884SHong Zhang 
1065801fbe65SHong Zhang     if (!mumps->myid) {
1066334c5f61SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,nrhs*M,&b_seq);CHKERRQ(ierr);
1067801fbe65SHong Zhang       ierr = ISCreateGeneral(PETSC_COMM_SELF,nrhs*M,idx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr);
1068801fbe65SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,nrhs*M,0,1,&is_from);CHKERRQ(ierr);
1069801fbe65SHong Zhang     } else {
1070334c5f61SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&b_seq);CHKERRQ(ierr);
1071801fbe65SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_to);CHKERRQ(ierr);
1072801fbe65SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_from);CHKERRQ(ierr);
1073801fbe65SHong Zhang     }
1074334c5f61SHong Zhang     ierr = VecScatterCreate(v_mpi,is_from,b_seq,is_to,&scat_rhs);CHKERRQ(ierr);
1075334c5f61SHong Zhang     ierr = VecScatterBegin(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1076801fbe65SHong Zhang     ierr = ISDestroy(&is_to);CHKERRQ(ierr);
1077801fbe65SHong Zhang     ierr = ISDestroy(&is_from);CHKERRQ(ierr);
1078334c5f61SHong Zhang     ierr = VecScatterEnd(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1079801fbe65SHong Zhang 
1080801fbe65SHong Zhang     if (!mumps->myid) { /* define rhs on the host */
1081334c5f61SHong Zhang       ierr = VecGetArray(b_seq,&bray);CHKERRQ(ierr);
1082940cd9d6SSatish Balay       mumps->id.rhs = (MumpsScalar*)bray;
1083334c5f61SHong Zhang       ierr = VecRestoreArray(b_seq,&bray);CHKERRQ(ierr);
1084801fbe65SHong Zhang     }
1085801fbe65SHong Zhang 
1086801fbe65SHong Zhang     /* solve phase */
1087801fbe65SHong Zhang     /*-------------*/
1088801fbe65SHong Zhang     mumps->id.job = JOB_SOLVE;
1089801fbe65SHong Zhang     PetscMUMPS_c(&mumps->id);
1090801fbe65SHong Zhang     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));
1091801fbe65SHong Zhang 
1092334c5f61SHong Zhang     /* scatter mumps distributed solution to petsc vector v_mpi, which shares local arrays with solution matrix X */
109374f0fcc7SHong Zhang     ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr);
109474f0fcc7SHong Zhang     ierr = VecPlaceArray(v_mpi,array);CHKERRQ(ierr);
1095801fbe65SHong Zhang 
1096334c5f61SHong Zhang     /* create scatter scat_sol */
109771aed81dSHong Zhang     ierr = PetscMalloc1(nlsol_loc,&idxx);CHKERRQ(ierr);
109871aed81dSHong Zhang     ierr = ISCreateStride(PETSC_COMM_SELF,nlsol_loc,0,1,&is_from);CHKERRQ(ierr);
109971aed81dSHong Zhang     for (i=0; i<lsol_loc; i++) {
1100334c5f61SHong Zhang       isol_loc[i] -= 1; /* change Fortran style to C style */
1101334c5f61SHong Zhang       idxx[i] = iidx[isol_loc[i]];
1102801fbe65SHong Zhang       for (j=1; j<nrhs; j++){
1103334c5f61SHong Zhang         idxx[j*lsol_loc+i] = iidx[isol_loc[i]+j*M];
1104801fbe65SHong Zhang       }
1105801fbe65SHong Zhang     }
110671aed81dSHong Zhang     ierr = ISCreateGeneral(PETSC_COMM_SELF,nlsol_loc,idxx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr);
1107334c5f61SHong Zhang     ierr = VecScatterCreate(x_seq,is_from,v_mpi,is_to,&scat_sol);CHKERRQ(ierr);
1108334c5f61SHong Zhang     ierr = VecScatterBegin(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1109801fbe65SHong Zhang     ierr = ISDestroy(&is_from);CHKERRQ(ierr);
1110801fbe65SHong Zhang     ierr = ISDestroy(&is_to);CHKERRQ(ierr);
1111334c5f61SHong Zhang     ierr = VecScatterEnd(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1112801fbe65SHong Zhang     ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr);
111371aed81dSHong Zhang 
111471aed81dSHong Zhang     /* free spaces */
111571aed81dSHong Zhang     mumps->id.sol_loc = sol_loc_save;
111671aed81dSHong Zhang     mumps->id.isol_loc = isol_loc_save;
111771aed81dSHong Zhang 
111871aed81dSHong Zhang     ierr = PetscFree2(sol_loc,isol_loc);CHKERRQ(ierr);
1119801fbe65SHong Zhang     ierr = PetscFree2(idx,iidx);CHKERRQ(ierr);
1120801fbe65SHong Zhang     ierr = PetscFree(idxx);CHKERRQ(ierr);
112171aed81dSHong Zhang     ierr = VecDestroy(&x_seq);CHKERRQ(ierr);
112274f0fcc7SHong Zhang     ierr = VecDestroy(&v_mpi);CHKERRQ(ierr);
1123334c5f61SHong Zhang     ierr = VecDestroy(&b_seq);CHKERRQ(ierr);
1124334c5f61SHong Zhang     ierr = VecScatterDestroy(&scat_rhs);CHKERRQ(ierr);
1125334c5f61SHong Zhang     ierr = VecScatterDestroy(&scat_sol);CHKERRQ(ierr);
1126801fbe65SHong Zhang   }
1127e0b74bf9SHong Zhang   PetscFunctionReturn(0);
1128e0b74bf9SHong Zhang }
1129e0b74bf9SHong Zhang 
1130ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX)
1131a58c3f20SHong Zhang /*
1132a58c3f20SHong Zhang   input:
1133a58c3f20SHong Zhang    F:        numeric factor
1134a58c3f20SHong Zhang   output:
1135a58c3f20SHong Zhang    nneg:     total number of negative pivots
113619d49a3bSHong Zhang    nzero:    total number of zero pivots
113719d49a3bSHong Zhang    npos:     (global dimension of F) - nneg - nzero
1138a58c3f20SHong Zhang */
1139dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
1140a58c3f20SHong Zhang {
1141e69c285eSBarry Smith   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1142dfbe8321SBarry Smith   PetscErrorCode ierr;
1143c1490034SHong Zhang   PetscMPIInt    size;
1144a58c3f20SHong Zhang 
1145a58c3f20SHong Zhang   PetscFunctionBegin;
1146ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)F),&size);CHKERRQ(ierr);
1147bcb30aebSHong Zhang   /* 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 */
1148a5e57a09SHong Zhang   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));
1149ed85ac9fSHong Zhang 
1150710ac8efSHong Zhang   if (nneg) *nneg = mumps->id.INFOG(12);
1151ed85ac9fSHong Zhang   if (nzero || npos) {
1152ed85ac9fSHong Zhang     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");
1153710ac8efSHong Zhang     if (nzero) *nzero = mumps->id.INFOG(28);
1154710ac8efSHong Zhang     if (npos) *npos   = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28));
1155a58c3f20SHong Zhang   }
1156a58c3f20SHong Zhang   PetscFunctionReturn(0);
1157a58c3f20SHong Zhang }
115819d49a3bSHong Zhang #endif
1159a58c3f20SHong Zhang 
11600481f469SBarry Smith PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
1161af281ebdSHong Zhang {
1162e69c285eSBarry Smith   Mat_MUMPS      *mumps =(Mat_MUMPS*)(F)->data;
11636849ba73SBarry Smith   PetscErrorCode ierr;
1164ace3abfcSBarry Smith   PetscBool      isMPIAIJ;
1165397b6df1SKris Buschelman 
1166397b6df1SKris Buschelman   PetscFunctionBegin;
11676baea169SHong Zhang   if (mumps->id.INFOG(1) < 0) {
11682aca8efcSHong Zhang     if (mumps->id.INFOG(1) == -6) {
11692aca8efcSHong Zhang       ierr = PetscInfo2(A,"MatFactorNumeric is called with singular matrix structure, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
11706baea169SHong Zhang     }
11716baea169SHong Zhang     ierr = PetscInfo2(A,"MatFactorNumeric is called after analysis phase fails, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
11722aca8efcSHong Zhang     PetscFunctionReturn(0);
11732aca8efcSHong Zhang   }
11746baea169SHong Zhang 
1175a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1176397b6df1SKris Buschelman 
1177397b6df1SKris Buschelman   /* numerical factorization phase */
1178329ec9b3SHong Zhang   /*-------------------------------*/
1179a5e57a09SHong Zhang   mumps->id.job = JOB_FACTNUMERIC;
11804e34a73bSHong Zhang   if (!mumps->id.ICNTL(18)) { /* A is centralized */
1181a5e57a09SHong Zhang     if (!mumps->myid) {
1182940cd9d6SSatish Balay       mumps->id.a = (MumpsScalar*)mumps->val;
1183397b6df1SKris Buschelman     }
1184397b6df1SKris Buschelman   } else {
1185940cd9d6SSatish Balay     mumps->id.a_loc = (MumpsScalar*)mumps->val;
1186397b6df1SKris Buschelman   }
1187a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
1188a5e57a09SHong Zhang   if (mumps->id.INFOG(1) < 0) {
1189c0d63f2fSHong Zhang     if (A->erroriffailure) {
1190c0d63f2fSHong Zhang       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1191151787a6SHong Zhang     } else {
1192c0d63f2fSHong Zhang       if (mumps->id.INFOG(1) == -10) { /* numerically singular matrix */
11932aca8efcSHong Zhang         ierr = PetscInfo2(F,"matrix is numerically singular, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
1194603e8f96SBarry Smith         F->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1195c0d63f2fSHong Zhang       } else if (mumps->id.INFOG(1) == -13) {
1196c0d63f2fSHong Zhang         ierr = PetscInfo2(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, cannot allocate required memory %d megabytes\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
1197603e8f96SBarry Smith         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1198c0d63f2fSHong Zhang       } else if (mumps->id.INFOG(1) == -8 || mumps->id.INFOG(1) == -9 || (-16 < mumps->id.INFOG(1) && mumps->id.INFOG(1) < -10) ) {
1199c0d63f2fSHong Zhang         ierr = PetscInfo2(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d, problem with workarray \n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
1200603e8f96SBarry Smith         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
12012aca8efcSHong Zhang       } else {
1202c0d63f2fSHong Zhang         ierr = PetscInfo2(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
1203603e8f96SBarry Smith         F->factorerrortype = MAT_FACTOR_OTHER;
1204151787a6SHong Zhang       }
12052aca8efcSHong Zhang     }
1206397b6df1SKris Buschelman   }
1207a5e57a09SHong Zhang   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));
1208397b6df1SKris Buschelman 
1209dcd589f8SShri Abhyankar   (F)->assembled        = PETSC_TRUE;
1210a5e57a09SHong Zhang   mumps->matstruc       = SAME_NONZERO_PATTERN;
1211b5fa320bSStefano Zampini   mumps->schur_factored = PETSC_FALSE;
121259ac8732SStefano Zampini   mumps->schur_inverted = PETSC_FALSE;
121367877ebaSShri Abhyankar 
1214066565c5SStefano Zampini   /* just to be sure that ICNTL(19) value returned by a call from MatMumpsGetIcntl is always consistent */
1215066565c5SStefano Zampini   if (!mumps->sym && mumps->id.ICNTL(19) && mumps->id.ICNTL(19) != 1) mumps->id.ICNTL(19) = 3;
1216066565c5SStefano Zampini 
1217a5e57a09SHong Zhang   if (mumps->size > 1) {
121867877ebaSShri Abhyankar     PetscInt    lsol_loc;
121967877ebaSShri Abhyankar     PetscScalar *sol_loc;
12202205254eSKarl Rupp 
1221c2093ab7SHong Zhang     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr);
1222c2093ab7SHong Zhang 
1223c2093ab7SHong Zhang     /* distributed solution; Create x_seq=sol_loc for repeated use */
1224c2093ab7SHong Zhang     if (mumps->x_seq) {
1225c2093ab7SHong Zhang       ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
1226c2093ab7SHong Zhang       ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr);
1227c2093ab7SHong Zhang       ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr);
1228c2093ab7SHong Zhang     }
1229a5e57a09SHong Zhang     lsol_loc = mumps->id.INFO(23); /* length of sol_loc */
1230dcca6d9dSJed Brown     ierr = PetscMalloc2(lsol_loc,&sol_loc,lsol_loc,&mumps->id.isol_loc);CHKERRQ(ierr);
1231a5e57a09SHong Zhang     mumps->id.lsol_loc = lsol_loc;
1232940cd9d6SSatish Balay     mumps->id.sol_loc = (MumpsScalar*)sol_loc;
1233a5e57a09SHong Zhang     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq);CHKERRQ(ierr);
123467877ebaSShri Abhyankar   }
1235397b6df1SKris Buschelman   PetscFunctionReturn(0);
1236397b6df1SKris Buschelman }
1237397b6df1SKris Buschelman 
12389a2535b5SHong Zhang /* Sets MUMPS options from the options database */
12399a2535b5SHong Zhang PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A)
1240dcd589f8SShri Abhyankar {
1241e69c285eSBarry Smith   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1242dcd589f8SShri Abhyankar   PetscErrorCode ierr;
1243b34f08ffSHong Zhang   PetscInt       icntl,info[40],i,ninfo=40;
1244ace3abfcSBarry Smith   PetscBool      flg;
1245dcd589f8SShri Abhyankar 
1246dcd589f8SShri Abhyankar   PetscFunctionBegin;
1247ce94432eSBarry Smith   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
12489a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);CHKERRQ(ierr);
12499a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(1) = icntl;
12509a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_2","ICNTL(2): output stream for diagnostic printing, statistics, and warning","None",mumps->id.ICNTL(2),&icntl,&flg);CHKERRQ(ierr);
12519a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(2) = icntl;
12529a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_3","ICNTL(3): output stream for global information, collected on the host","None",mumps->id.ICNTL(3),&icntl,&flg);CHKERRQ(ierr);
12539a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(3) = icntl;
1254dcd589f8SShri Abhyankar 
12559a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",mumps->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
12569a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(4) = icntl;
12579a2535b5SHong Zhang   if (mumps->id.ICNTL(4) || PetscLogPrintInfo) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */
12589a2535b5SHong Zhang 
1259d341cd04SHong Zhang   ierr = 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);CHKERRQ(ierr);
12609a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(6) = icntl;
12619a2535b5SHong Zhang 
1262d341cd04SHong Zhang   ierr = 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);CHKERRQ(ierr);
1263dcd589f8SShri Abhyankar   if (flg) {
12642205254eSKarl Rupp     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");
12652205254eSKarl Rupp     else mumps->id.ICNTL(7) = icntl;
1266dcd589f8SShri Abhyankar   }
1267e0b74bf9SHong Zhang 
12680298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_8","ICNTL(8): scaling strategy (-2 to 8 or 77)","None",mumps->id.ICNTL(8),&mumps->id.ICNTL(8),NULL);CHKERRQ(ierr);
1269d341cd04SHong Zhang   /* ierr = 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);CHKERRQ(ierr); handled by MatSolveTranspose_MUMPS() */
12700298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",mumps->id.ICNTL(10),&mumps->id.ICNTL(10),NULL);CHKERRQ(ierr);
1271d341cd04SHong Zhang   ierr = 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);CHKERRQ(ierr);
1272d341cd04SHong Zhang   ierr = 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);CHKERRQ(ierr);
1273d341cd04SHong Zhang   ierr = 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);CHKERRQ(ierr);
1274d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage increase in the estimated working space","None",mumps->id.ICNTL(14),&mumps->id.ICNTL(14),NULL);CHKERRQ(ierr);
1275d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): computes the Schur complement","None",mumps->id.ICNTL(19),&mumps->id.ICNTL(19),NULL);CHKERRQ(ierr);
127659ac8732SStefano Zampini   if (mumps->id.ICNTL(19) <= 0 || mumps->id.ICNTL(19) > 3) { /* reset any schur data (if any) */
127759ac8732SStefano Zampini     ierr = MatMumpsResetSchur_Private(mumps);CHKERRQ(ierr);
127859ac8732SStefano Zampini   }
12794e34a73bSHong Zhang   /* ierr = 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);CHKERRQ(ierr); -- sparse rhs is not supported in PETSc API */
1280d341cd04SHong Zhang   /* ierr = 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);CHKERRQ(ierr); we only use distributed solution vector */
12819a2535b5SHong Zhang 
1282d341cd04SHong Zhang   ierr = 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);CHKERRQ(ierr);
12830298fd71SBarry Smith   ierr = 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);CHKERRQ(ierr);
12840298fd71SBarry Smith   ierr = 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);CHKERRQ(ierr);
12859a2535b5SHong Zhang   if (mumps->id.ICNTL(24)) {
12869a2535b5SHong Zhang     mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */
1287d7ebd59bSHong Zhang   }
1288d7ebd59bSHong Zhang 
1289d341cd04SHong Zhang   ierr = 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);CHKERRQ(ierr);
1290d341cd04SHong Zhang   ierr = 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);CHKERRQ(ierr);
12912cd7d884SHong Zhang   ierr = 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);CHKERRQ(ierr);
12920298fd71SBarry Smith   ierr = 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);CHKERRQ(ierr);
1293d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_29","ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis","None",mumps->id.ICNTL(29),&mumps->id.ICNTL(29),NULL);CHKERRQ(ierr);
12940298fd71SBarry Smith   ierr = 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);CHKERRQ(ierr);
1295d341cd04SHong Zhang   ierr = 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);CHKERRQ(ierr);
12964e34a73bSHong Zhang   /* ierr = 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);CHKERRQ(ierr);  -- not supported by PETSc API */
12970298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),NULL);CHKERRQ(ierr);
1298dcd589f8SShri Abhyankar 
12990298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),NULL);CHKERRQ(ierr);
13000298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",mumps->id.CNTL(2),&mumps->id.CNTL(2),NULL);CHKERRQ(ierr);
13010298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),NULL);CHKERRQ(ierr);
13020298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",mumps->id.CNTL(4),&mumps->id.CNTL(4),NULL);CHKERRQ(ierr);
13030298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",mumps->id.CNTL(5),&mumps->id.CNTL(5),NULL);CHKERRQ(ierr);
1304e5bb22a1SHong Zhang 
13052a808120SBarry Smith   ierr = PetscOptionsString("-mat_mumps_ooc_tmpdir", "out of core directory", "None", mumps->id.ooc_tmpdir, mumps->id.ooc_tmpdir, 256, NULL);CHKERRQ(ierr);
1306b34f08ffSHong Zhang 
130716d797efSHong Zhang   ierr = PetscOptionsIntArray("-mat_mumps_view_info","request INFO local to each processor","",info,&ninfo,NULL);CHKERRQ(ierr);
1308b34f08ffSHong Zhang   if (ninfo) {
1309b34f08ffSHong Zhang     if (ninfo > 40) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"number of INFO %d must <= 40\n",ninfo);
1310b34f08ffSHong Zhang     ierr = PetscMalloc1(ninfo,&mumps->info);CHKERRQ(ierr);
1311b34f08ffSHong Zhang     mumps->ninfo = ninfo;
1312b34f08ffSHong Zhang     for (i=0; i<ninfo; i++) {
13136c4ed002SBarry Smith       if (info[i] < 0 || info[i]>40) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index of INFO %d must between 1 and 40\n",ninfo);
13142a808120SBarry Smith       else  mumps->info[i] = info[i];
1315b34f08ffSHong Zhang     }
1316b34f08ffSHong Zhang   }
1317b34f08ffSHong Zhang 
13182a808120SBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1319dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
1320dcd589f8SShri Abhyankar }
1321dcd589f8SShri Abhyankar 
1322f697e70eSHong Zhang PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps)
1323dcd589f8SShri Abhyankar {
1324dcd589f8SShri Abhyankar   PetscErrorCode ierr;
1325dcd589f8SShri Abhyankar 
1326dcd589f8SShri Abhyankar   PetscFunctionBegin;
13272a808120SBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)A), &mumps->myid);CHKERRQ(ierr);
1328ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&mumps->size);CHKERRQ(ierr);
1329ce94432eSBarry Smith   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(mumps->comm_mumps));CHKERRQ(ierr);
13302205254eSKarl Rupp 
1331f697e70eSHong Zhang   mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps);
1332f697e70eSHong Zhang 
1333f697e70eSHong Zhang   mumps->id.job = JOB_INIT;
1334f697e70eSHong Zhang   mumps->id.par = 1;  /* host participates factorizaton and solve */
1335f697e70eSHong Zhang   mumps->id.sym = mumps->sym;
13362907cef9SHong Zhang   PetscMUMPS_c(&mumps->id);
1337f697e70eSHong Zhang 
13380298fd71SBarry Smith   mumps->scat_rhs     = NULL;
13390298fd71SBarry Smith   mumps->scat_sol     = NULL;
13409a2535b5SHong Zhang 
134170544d5fSHong Zhang   /* set PETSc-MUMPS default options - override MUMPS default */
13429a2535b5SHong Zhang   mumps->id.ICNTL(3) = 0;
13439a2535b5SHong Zhang   mumps->id.ICNTL(4) = 0;
13449a2535b5SHong Zhang   if (mumps->size == 1) {
13459a2535b5SHong Zhang     mumps->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
13469a2535b5SHong Zhang   } else {
13479a2535b5SHong Zhang     mumps->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
13484e34a73bSHong Zhang     mumps->id.ICNTL(20) = 0;   /* rhs is in dense format */
134970544d5fSHong Zhang     mumps->id.ICNTL(21) = 1;   /* distributed solution */
13509a2535b5SHong Zhang   }
13516444a565SStefano Zampini 
13526444a565SStefano Zampini   /* schur */
13536444a565SStefano Zampini   mumps->id.size_schur      = 0;
13546444a565SStefano Zampini   mumps->id.listvar_schur   = NULL;
13556444a565SStefano Zampini   mumps->id.schur           = NULL;
1356b5fa320bSStefano Zampini   mumps->sizeredrhs         = 0;
1357b5fa320bSStefano Zampini   mumps->schur_pivots       = NULL;
1358b5fa320bSStefano Zampini   mumps->schur_work         = NULL;
135959ac8732SStefano Zampini   mumps->schur_sol          = NULL;
136059ac8732SStefano Zampini   mumps->schur_sizesol      = 0;
136159ac8732SStefano Zampini   mumps->schur_factored     = PETSC_FALSE;
136259ac8732SStefano Zampini   mumps->schur_inverted     = PETSC_FALSE;
1363a0b0af32SStefano Zampini   mumps->schur_sym          = mumps->id.sym;
1364dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
1365dcd589f8SShri Abhyankar }
1366dcd589f8SShri Abhyankar 
13679a625307SHong Zhang PetscErrorCode MatFactorSymbolic_MUMPS_ReportIfError(Mat F,Mat A,const MatFactorInfo *info,Mat_MUMPS *mumps)
13685cd7cf9dSHong Zhang {
13695cd7cf9dSHong Zhang   PetscErrorCode ierr;
13705cd7cf9dSHong Zhang 
13715cd7cf9dSHong Zhang   PetscFunctionBegin;
13725cd7cf9dSHong Zhang   if (mumps->id.INFOG(1) < 0) {
13735cd7cf9dSHong Zhang     if (A->erroriffailure) {
13745cd7cf9dSHong Zhang       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
13755cd7cf9dSHong Zhang     } else {
13765cd7cf9dSHong Zhang       if (mumps->id.INFOG(1) == -6) {
13775cd7cf9dSHong Zhang         ierr = PetscInfo2(F,"matrix is singular in structure, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
1378603e8f96SBarry Smith         F->factorerrortype = MAT_FACTOR_STRUCT_ZEROPIVOT;
13795cd7cf9dSHong Zhang       } else if (mumps->id.INFOG(1) == -5 || mumps->id.INFOG(1) == -7) {
13805cd7cf9dSHong Zhang         ierr = PetscInfo2(F,"problem of workspace, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
1381603e8f96SBarry Smith         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
13825cd7cf9dSHong Zhang       } else {
13835cd7cf9dSHong Zhang         ierr = PetscInfo2(F,"Error reported by MUMPS in analysis phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
1384603e8f96SBarry Smith         F->factorerrortype = MAT_FACTOR_OTHER;
13855cd7cf9dSHong Zhang       }
13865cd7cf9dSHong Zhang     }
13875cd7cf9dSHong Zhang   }
13885cd7cf9dSHong Zhang   PetscFunctionReturn(0);
13895cd7cf9dSHong Zhang }
13905cd7cf9dSHong Zhang 
1391a5e57a09SHong Zhang /* Note Petsc r(=c) permutation is used when mumps->id.ICNTL(7)==1 with centralized assembled matrix input; otherwise r and c are ignored */
13920481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1393b24902e0SBarry Smith {
1394e69c285eSBarry Smith   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1395dcd589f8SShri Abhyankar   PetscErrorCode ierr;
139667877ebaSShri Abhyankar   Vec            b;
139767877ebaSShri Abhyankar   IS             is_iden;
139867877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1399397b6df1SKris Buschelman 
1400397b6df1SKris Buschelman   PetscFunctionBegin;
1401a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1402dcd589f8SShri Abhyankar 
14039a2535b5SHong Zhang   /* Set MUMPS options from the options database */
14049a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1405dcd589f8SShri Abhyankar 
1406a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1407dcd589f8SShri Abhyankar 
140867877ebaSShri Abhyankar   /* analysis phase */
140967877ebaSShri Abhyankar   /*----------------*/
1410a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1411a5e57a09SHong Zhang   mumps->id.n   = M;
1412a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
141367877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1414a5e57a09SHong Zhang     if (!mumps->myid) {
1415a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1416a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
1417940cd9d6SSatish Balay         mumps->id.a = (MumpsScalar*)mumps->val;
141867877ebaSShri Abhyankar       }
1419a5e57a09SHong Zhang       if (mumps->id.ICNTL(7) == 1) { /* use user-provide matrix ordering - assuming r = c ordering */
14205248a706SHong Zhang         /*
14215248a706SHong Zhang         PetscBool      flag;
14225248a706SHong Zhang         ierr = ISEqual(r,c,&flag);CHKERRQ(ierr);
14235248a706SHong Zhang         if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"row_perm != col_perm");
14245248a706SHong Zhang         ierr = ISView(r,PETSC_VIEWER_STDOUT_SELF);
14255248a706SHong Zhang          */
1426a5e57a09SHong Zhang         if (!mumps->myid) {
1427e0b74bf9SHong Zhang           const PetscInt *idx;
1428e0b74bf9SHong Zhang           PetscInt       i,*perm_in;
14292205254eSKarl Rupp 
1430785e854fSJed Brown           ierr = PetscMalloc1(M,&perm_in);CHKERRQ(ierr);
1431e0b74bf9SHong Zhang           ierr = ISGetIndices(r,&idx);CHKERRQ(ierr);
14322205254eSKarl Rupp 
1433a5e57a09SHong Zhang           mumps->id.perm_in = perm_in;
1434e0b74bf9SHong Zhang           for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */
1435e0b74bf9SHong Zhang           ierr = ISRestoreIndices(r,&idx);CHKERRQ(ierr);
1436e0b74bf9SHong Zhang         }
1437e0b74bf9SHong Zhang       }
143867877ebaSShri Abhyankar     }
143967877ebaSShri Abhyankar     break;
144067877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1441a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1442a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1443a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
1444940cd9d6SSatish Balay       mumps->id.a_loc = (MumpsScalar*)mumps->val;
144567877ebaSShri Abhyankar     }
144667877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1447a5e57a09SHong Zhang     if (!mumps->myid) {
14482cd7d884SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->rmap->N,&mumps->b_seq);CHKERRQ(ierr);
14492cd7d884SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->N,0,1,&is_iden);CHKERRQ(ierr);
145067877ebaSShri Abhyankar     } else {
1451a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
145267877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
145367877ebaSShri Abhyankar     }
14542a7a6963SBarry Smith     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1455a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
14566bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
14576bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
145867877ebaSShri Abhyankar     break;
145967877ebaSShri Abhyankar   }
1460a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
14615cd7cf9dSHong Zhang   ierr = MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);CHKERRQ(ierr);
146267877ebaSShri Abhyankar 
1463719d5645SBarry Smith   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1464dcd589f8SShri Abhyankar   F->ops->solve           = MatSolve_MUMPS;
146551d5961aSHong Zhang   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
14664e34a73bSHong Zhang   F->ops->matsolve        = MatMatSolve_MUMPS;
1467b24902e0SBarry Smith   PetscFunctionReturn(0);
1468b24902e0SBarry Smith }
1469b24902e0SBarry Smith 
1470450b117fSShri Abhyankar /* Note the Petsc r and c permutations are ignored */
1471450b117fSShri Abhyankar PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1472450b117fSShri Abhyankar {
1473e69c285eSBarry Smith   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1474dcd589f8SShri Abhyankar   PetscErrorCode ierr;
147567877ebaSShri Abhyankar   Vec            b;
147667877ebaSShri Abhyankar   IS             is_iden;
147767877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1478450b117fSShri Abhyankar 
1479450b117fSShri Abhyankar   PetscFunctionBegin;
1480a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1481dcd589f8SShri Abhyankar 
14829a2535b5SHong Zhang   /* Set MUMPS options from the options database */
14839a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1484dcd589f8SShri Abhyankar 
1485a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
148667877ebaSShri Abhyankar 
148767877ebaSShri Abhyankar   /* analysis phase */
148867877ebaSShri Abhyankar   /*----------------*/
1489a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1490a5e57a09SHong Zhang   mumps->id.n   = M;
1491a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
149267877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1493a5e57a09SHong Zhang     if (!mumps->myid) {
1494a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1495a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
1496940cd9d6SSatish Balay         mumps->id.a = (MumpsScalar*)mumps->val;
149767877ebaSShri Abhyankar       }
149867877ebaSShri Abhyankar     }
149967877ebaSShri Abhyankar     break;
150067877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1501a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1502a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1503a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
1504940cd9d6SSatish Balay       mumps->id.a_loc = (MumpsScalar*)mumps->val;
150567877ebaSShri Abhyankar     }
150667877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1507a5e57a09SHong Zhang     if (!mumps->myid) {
1508a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
150967877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
151067877ebaSShri Abhyankar     } else {
1511a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
151267877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
151367877ebaSShri Abhyankar     }
15142a7a6963SBarry Smith     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1515a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
15166bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
15176bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
151867877ebaSShri Abhyankar     break;
151967877ebaSShri Abhyankar   }
1520a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
15215cd7cf9dSHong Zhang   ierr = MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);CHKERRQ(ierr);
152267877ebaSShri Abhyankar 
1523450b117fSShri Abhyankar   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1524dcd589f8SShri Abhyankar   F->ops->solve           = MatSolve_MUMPS;
152551d5961aSHong Zhang   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1526450b117fSShri Abhyankar   PetscFunctionReturn(0);
1527450b117fSShri Abhyankar }
1528b24902e0SBarry Smith 
1529141f4205SHong Zhang /* Note the Petsc r permutation and factor info are ignored */
153067877ebaSShri Abhyankar PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
1531b24902e0SBarry Smith {
1532e69c285eSBarry Smith   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1533dcd589f8SShri Abhyankar   PetscErrorCode ierr;
153467877ebaSShri Abhyankar   Vec            b;
153567877ebaSShri Abhyankar   IS             is_iden;
153667877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1537397b6df1SKris Buschelman 
1538397b6df1SKris Buschelman   PetscFunctionBegin;
1539a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1540dcd589f8SShri Abhyankar 
15419a2535b5SHong Zhang   /* Set MUMPS options from the options database */
15429a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1543dcd589f8SShri Abhyankar 
1544a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1545dcd589f8SShri Abhyankar 
154667877ebaSShri Abhyankar   /* analysis phase */
154767877ebaSShri Abhyankar   /*----------------*/
1548a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1549a5e57a09SHong Zhang   mumps->id.n   = M;
1550a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
155167877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1552a5e57a09SHong Zhang     if (!mumps->myid) {
1553a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1554a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
1555940cd9d6SSatish Balay         mumps->id.a = (MumpsScalar*)mumps->val;
155667877ebaSShri Abhyankar       }
155767877ebaSShri Abhyankar     }
155867877ebaSShri Abhyankar     break;
155967877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1560a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1561a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1562a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
1563940cd9d6SSatish Balay       mumps->id.a_loc = (MumpsScalar*)mumps->val;
156467877ebaSShri Abhyankar     }
156567877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1566a5e57a09SHong Zhang     if (!mumps->myid) {
1567a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
156867877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
156967877ebaSShri Abhyankar     } else {
1570a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
157167877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
157267877ebaSShri Abhyankar     }
15732a7a6963SBarry Smith     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1574a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
15756bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
15766bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
157767877ebaSShri Abhyankar     break;
157867877ebaSShri Abhyankar   }
1579a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
15805cd7cf9dSHong Zhang   ierr = MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);CHKERRQ(ierr);
15815cd7cf9dSHong Zhang 
15822792810eSHong Zhang   F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
1583dcd589f8SShri Abhyankar   F->ops->solve                 = MatSolve_MUMPS;
158451d5961aSHong Zhang   F->ops->solvetranspose        = MatSolve_MUMPS;
15854e34a73bSHong Zhang   F->ops->matsolve              = MatMatSolve_MUMPS;
15864e34a73bSHong Zhang #if defined(PETSC_USE_COMPLEX)
15870298fd71SBarry Smith   F->ops->getinertia = NULL;
15884e34a73bSHong Zhang #else
15894e34a73bSHong Zhang   F->ops->getinertia = MatGetInertia_SBAIJMUMPS;
1590db4efbfdSBarry Smith #endif
1591b24902e0SBarry Smith   PetscFunctionReturn(0);
1592b24902e0SBarry Smith }
1593b24902e0SBarry Smith 
159464e6c443SBarry Smith PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
159574ed9c26SBarry Smith {
1596f6c57405SHong Zhang   PetscErrorCode    ierr;
159764e6c443SBarry Smith   PetscBool         iascii;
159864e6c443SBarry Smith   PetscViewerFormat format;
1599e69c285eSBarry Smith   Mat_MUMPS         *mumps=(Mat_MUMPS*)A->data;
1600f6c57405SHong Zhang 
1601f6c57405SHong Zhang   PetscFunctionBegin;
160264e6c443SBarry Smith   /* check if matrix is mumps type */
160364e6c443SBarry Smith   if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0);
160464e6c443SBarry Smith 
1605251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
160664e6c443SBarry Smith   if (iascii) {
160764e6c443SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
160864e6c443SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO) {
160964e6c443SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
1610a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                   %d \n",mumps->id.sym);CHKERRQ(ierr);
1611a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):            %d \n",mumps->id.par);CHKERRQ(ierr);
1612a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):         %d \n",mumps->id.ICNTL(1));CHKERRQ(ierr);
1613a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg): %d \n",mumps->id.ICNTL(2));CHKERRQ(ierr);
1614a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):   %d \n",mumps->id.ICNTL(3));CHKERRQ(ierr);
1615a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):        %d \n",mumps->id.ICNTL(4));CHKERRQ(ierr);
1616a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):         %d \n",mumps->id.ICNTL(5));CHKERRQ(ierr);
1617a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):        %d \n",mumps->id.ICNTL(6));CHKERRQ(ierr);
1618a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (sequentia matrix ordering):%d \n",mumps->id.ICNTL(7));CHKERRQ(ierr);
1619a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scalling strategy):        %d \n",mumps->id.ICNTL(8));CHKERRQ(ierr);
1620a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements):  %d \n",mumps->id.ICNTL(10));CHKERRQ(ierr);
1621a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):          %d \n",mumps->id.ICNTL(11));CHKERRQ(ierr);
1622a5e57a09SHong Zhang       if (mumps->id.ICNTL(11)>0) {
1623a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(4) (inf norm of input mat):        %g\n",mumps->id.RINFOG(4));CHKERRQ(ierr);
1624a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(5) (inf norm of solution):         %g\n",mumps->id.RINFOG(5));CHKERRQ(ierr);
1625a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(6) (inf norm of residual):         %g\n",mumps->id.RINFOG(6));CHKERRQ(ierr);
1626a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8));CHKERRQ(ierr);
1627a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(9) (error estimate):               %g \n",mumps->id.RINFOG(9));CHKERRQ(ierr);
1628a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11));CHKERRQ(ierr);
1629f6c57405SHong Zhang       }
1630a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",mumps->id.ICNTL(12));CHKERRQ(ierr);
1631a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",mumps->id.ICNTL(13));CHKERRQ(ierr);
1632a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",mumps->id.ICNTL(14));CHKERRQ(ierr);
1633f6c57405SHong Zhang       /* ICNTL(15-17) not used */
1634a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",mumps->id.ICNTL(18));CHKERRQ(ierr);
1635a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Shur complement info):                       %d \n",mumps->id.ICNTL(19));CHKERRQ(ierr);
1636a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",mumps->id.ICNTL(20));CHKERRQ(ierr);
1637ca3dc52bSPierre Jolivet       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d \n",mumps->id.ICNTL(21));CHKERRQ(ierr);
1638a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",mumps->id.ICNTL(22));CHKERRQ(ierr);
1639a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",mumps->id.ICNTL(23));CHKERRQ(ierr);
1640c0165424SHong Zhang 
1641a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",mumps->id.ICNTL(24));CHKERRQ(ierr);
1642a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",mumps->id.ICNTL(25));CHKERRQ(ierr);
1643a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",mumps->id.ICNTL(26));CHKERRQ(ierr);
1644a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",mumps->id.ICNTL(27));CHKERRQ(ierr);
1645a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(28) (use parallel or sequential ordering):        %d \n",mumps->id.ICNTL(28));CHKERRQ(ierr);
1646a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(29) (parallel ordering):                          %d \n",mumps->id.ICNTL(29));CHKERRQ(ierr);
164742179a6aSHong Zhang 
1648a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(30) (user-specified set of entries in inv(A)):    %d \n",mumps->id.ICNTL(30));CHKERRQ(ierr);
1649a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(31) (factors is discarded in the solve phase):    %d \n",mumps->id.ICNTL(31));CHKERRQ(ierr);
1650a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(33) (compute determinant):                        %d \n",mumps->id.ICNTL(33));CHKERRQ(ierr);
1651f6c57405SHong Zhang 
1652a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",mumps->id.CNTL(1));CHKERRQ(ierr);
1653a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2));CHKERRQ(ierr);
1654ca3dc52bSPierre Jolivet       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",mumps->id.CNTL(3));CHKERRQ(ierr);
1655ca3dc52bSPierre Jolivet       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",mumps->id.CNTL(4));CHKERRQ(ierr);
1656a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",mumps->id.CNTL(5));CHKERRQ(ierr);
1657f6c57405SHong Zhang 
1658f6c57405SHong Zhang       /* infomation local to each processor */
165934ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);
16601575c14dSBarry Smith       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
1661a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %g \n",mumps->myid,mumps->id.RINFO(1));CHKERRQ(ierr);
16622a808120SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
166334ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);
1664a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(2));CHKERRQ(ierr);
16652a808120SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
166634ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);
1667a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(3));CHKERRQ(ierr);
16682a808120SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1669f6c57405SHong Zhang 
167034ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr);
1671a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] %d \n",mumps->myid,mumps->id.INFO(15));CHKERRQ(ierr);
16722a808120SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1673f6c57405SHong Zhang 
167434ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr);
1675a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(16));CHKERRQ(ierr);
16762a808120SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1677f6c57405SHong Zhang 
167834ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);
1679a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(23));CHKERRQ(ierr);
16802a808120SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1681b34f08ffSHong Zhang 
1682b34f08ffSHong Zhang       if (mumps->ninfo && mumps->ninfo <= 40){
1683b34f08ffSHong Zhang         PetscInt i;
1684b34f08ffSHong Zhang         for (i=0; i<mumps->ninfo; i++){
1685b34f08ffSHong Zhang           ierr = PetscViewerASCIIPrintf(viewer, "  INFO(%d): \n",mumps->info[i]);CHKERRQ(ierr);
1686b34f08ffSHong Zhang           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(mumps->info[i]));CHKERRQ(ierr);
16872a808120SBarry Smith           ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1688b34f08ffSHong Zhang         }
1689b34f08ffSHong Zhang       }
1690b34f08ffSHong Zhang 
1691b34f08ffSHong Zhang 
16921575c14dSBarry Smith       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
1693f6c57405SHong Zhang 
1694a5e57a09SHong Zhang       if (!mumps->myid) { /* information from the host */
1695a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",mumps->id.RINFOG(1));CHKERRQ(ierr);
1696a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",mumps->id.RINFOG(2));CHKERRQ(ierr);
1697a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",mumps->id.RINFOG(3));CHKERRQ(ierr);
1698a5e57a09SHong Zhang         ierr = 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));CHKERRQ(ierr);
1699f6c57405SHong Zhang 
1700a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(3));CHKERRQ(ierr);
1701a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(4));CHKERRQ(ierr);
1702a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",mumps->id.INFOG(5));CHKERRQ(ierr);
1703a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",mumps->id.INFOG(6));CHKERRQ(ierr);
1704a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively use after analysis): %d \n",mumps->id.INFOG(7));CHKERRQ(ierr);
1705a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",mumps->id.INFOG(8));CHKERRQ(ierr);
1706a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",mumps->id.INFOG(9));CHKERRQ(ierr);
1707a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",mumps->id.INFOG(10));CHKERRQ(ierr);
1708a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",mumps->id.INFOG(11));CHKERRQ(ierr);
1709a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",mumps->id.INFOG(12));CHKERRQ(ierr);
1710a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",mumps->id.INFOG(13));CHKERRQ(ierr);
1711a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",mumps->id.INFOG(14));CHKERRQ(ierr);
1712a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",mumps->id.INFOG(15));CHKERRQ(ierr);
1713a5e57a09SHong Zhang         ierr = 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));CHKERRQ(ierr);
1714a5e57a09SHong Zhang         ierr = 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));CHKERRQ(ierr);
1715a5e57a09SHong Zhang         ierr = 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));CHKERRQ(ierr);
1716a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(19) (size of all MUMPS internal data allocated during factorization: sum over all processors): %d \n",mumps->id.INFOG(19));CHKERRQ(ierr);
1717a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",mumps->id.INFOG(20));CHKERRQ(ierr);
1718a5e57a09SHong Zhang         ierr = 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));CHKERRQ(ierr);
1719a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(22) (size in MB of memory effectively used during factorization - sum over all processors): %d \n",mumps->id.INFOG(22));CHKERRQ(ierr);
1720a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",mumps->id.INFOG(23));CHKERRQ(ierr);
1721a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",mumps->id.INFOG(24));CHKERRQ(ierr);
1722a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",mumps->id.INFOG(25));CHKERRQ(ierr);
172340d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(28) (after factorization: number of null pivots encountered): %d\n",mumps->id.INFOG(28));CHKERRQ(ierr);
172440d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(29) (after factorization: effective number of entries in the factors (sum over all processors)): %d\n",mumps->id.INFOG(29));CHKERRQ(ierr);
172540d435e3SHong Zhang         ierr = 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));CHKERRQ(ierr);
172640d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(32) (after analysis: type of analysis done): %d\n",mumps->id.INFOG(32));CHKERRQ(ierr);
172740d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(33) (value used for ICNTL(8)): %d\n",mumps->id.INFOG(33));CHKERRQ(ierr);
172840d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(34) (exponent of the determinant if determinant is requested): %d\n",mumps->id.INFOG(34));CHKERRQ(ierr);
1729f6c57405SHong Zhang       }
1730f6c57405SHong Zhang     }
1731cb828f0fSHong Zhang   }
1732f6c57405SHong Zhang   PetscFunctionReturn(0);
1733f6c57405SHong Zhang }
1734f6c57405SHong Zhang 
173535bd34faSBarry Smith PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
173635bd34faSBarry Smith {
1737e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)A->data;
173835bd34faSBarry Smith 
173935bd34faSBarry Smith   PetscFunctionBegin;
174035bd34faSBarry Smith   info->block_size        = 1.0;
1741cb828f0fSHong Zhang   info->nz_allocated      = mumps->id.INFOG(20);
1742cb828f0fSHong Zhang   info->nz_used           = mumps->id.INFOG(20);
174335bd34faSBarry Smith   info->nz_unneeded       = 0.0;
174435bd34faSBarry Smith   info->assemblies        = 0.0;
174535bd34faSBarry Smith   info->mallocs           = 0.0;
174635bd34faSBarry Smith   info->memory            = 0.0;
174735bd34faSBarry Smith   info->fill_ratio_given  = 0;
174835bd34faSBarry Smith   info->fill_ratio_needed = 0;
174935bd34faSBarry Smith   info->factor_mallocs    = 0;
175035bd34faSBarry Smith   PetscFunctionReturn(0);
175135bd34faSBarry Smith }
175235bd34faSBarry Smith 
17535ccb76cbSHong Zhang /* -------------------------------------------------------------------------------------------*/
17548e7ba810SStefano Zampini PetscErrorCode MatFactorSetSchurIS_MUMPS(Mat F, IS is)
17556444a565SStefano Zampini {
1756e69c285eSBarry Smith   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
17578e7ba810SStefano Zampini   const PetscInt *idxs;
17588e7ba810SStefano Zampini   PetscInt       size,i;
17596444a565SStefano Zampini   PetscErrorCode ierr;
17606444a565SStefano Zampini 
17616444a565SStefano Zampini   PetscFunctionBegin;
1762241dbb5eSStefano Zampini   if (mumps->size > 1) {
1763241dbb5eSStefano Zampini     PetscBool ls,gs;
1764241dbb5eSStefano Zampini 
1765241dbb5eSStefano Zampini     ierr = ISGetLocalSize(is,&size);CHKERRQ(ierr);
17664c644ebcSSatish Balay     ls   = mumps->myid ? (size ? PETSC_FALSE : PETSC_TRUE) : PETSC_TRUE;
1767241dbb5eSStefano Zampini     ierr = MPI_Allreduce(&ls,&gs,1,MPIU_BOOL,MPI_LAND,mumps->comm_mumps);CHKERRQ(ierr);
1768241dbb5eSStefano Zampini     if (!gs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MUMPS distributed parallel Schur complements not yet supported from PETSc\n");
1769241dbb5eSStefano Zampini   }
17708e7ba810SStefano Zampini   ierr = ISGetLocalSize(is,&size);CHKERRQ(ierr);
17716444a565SStefano Zampini   if (mumps->id.size_schur != size) {
17726444a565SStefano Zampini     ierr = PetscFree2(mumps->id.listvar_schur,mumps->id.schur);CHKERRQ(ierr);
17736444a565SStefano Zampini     mumps->id.size_schur = size;
17746444a565SStefano Zampini     mumps->id.schur_lld = size;
17756444a565SStefano Zampini     ierr = PetscMalloc2(size,&mumps->id.listvar_schur,size*size,&mumps->id.schur);CHKERRQ(ierr);
17766444a565SStefano Zampini   }
17778e7ba810SStefano Zampini   ierr = ISGetIndices(is,&idxs);CHKERRQ(ierr);
17786444a565SStefano Zampini   ierr = PetscMemcpy(mumps->id.listvar_schur,idxs,size*sizeof(PetscInt));CHKERRQ(ierr);
17798e7ba810SStefano Zampini   /* MUMPS expects Fortran style indices */
17808e7ba810SStefano Zampini   for (i=0;i<size;i++) mumps->id.listvar_schur[i]++;
17818e7ba810SStefano Zampini   ierr = ISRestoreIndices(is,&idxs);CHKERRQ(ierr);
1782241dbb5eSStefano Zampini   if (mumps->size > 1) {
1783241dbb5eSStefano Zampini     mumps->id.ICNTL(19) = 1; /* MUMPS returns Schur centralized on the host */
1784241dbb5eSStefano Zampini   } else {
17856444a565SStefano Zampini     if (F->factortype == MAT_FACTOR_LU) {
178659ac8732SStefano Zampini       mumps->id.ICNTL(19) = 3; /* MUMPS returns full matrix */
17876444a565SStefano Zampini     } else {
178859ac8732SStefano Zampini       mumps->id.ICNTL(19) = 2; /* MUMPS returns lower triangular part */
17896444a565SStefano Zampini     }
1790241dbb5eSStefano Zampini   }
179159ac8732SStefano Zampini   /* set a special value of ICNTL (not handled my MUMPS) to be used in the solve phase by PETSc */
1792b5fa320bSStefano Zampini   mumps->id.ICNTL(26) = -1;
17936444a565SStefano Zampini   PetscFunctionReturn(0);
17946444a565SStefano Zampini }
179559ac8732SStefano Zampini 
17966444a565SStefano Zampini /* -------------------------------------------------------------------------------------------*/
17975a05ddb0SStefano Zampini PetscErrorCode MatFactorCreateSchurComplement_MUMPS(Mat F,Mat* S)
17986444a565SStefano Zampini {
17996444a565SStefano Zampini   Mat            St;
1800e69c285eSBarry Smith   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
18016444a565SStefano Zampini   PetscScalar    *array;
18026444a565SStefano Zampini #if defined(PETSC_USE_COMPLEX)
18038ac429a0SStefano Zampini   PetscScalar    im = PetscSqrtScalar((PetscScalar)-1.0);
18046444a565SStefano Zampini #endif
18056444a565SStefano Zampini   PetscErrorCode ierr;
18066444a565SStefano Zampini 
18076444a565SStefano Zampini   PetscFunctionBegin;
18085a05ddb0SStefano Zampini   if (!mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it");
1809241dbb5eSStefano Zampini   ierr = MatCreate(PETSC_COMM_SELF,&St);CHKERRQ(ierr);
18106444a565SStefano Zampini   ierr = MatSetSizes(St,PETSC_DECIDE,PETSC_DECIDE,mumps->id.size_schur,mumps->id.size_schur);CHKERRQ(ierr);
18116444a565SStefano Zampini   ierr = MatSetType(St,MATDENSE);CHKERRQ(ierr);
18126444a565SStefano Zampini   ierr = MatSetUp(St);CHKERRQ(ierr);
18136444a565SStefano Zampini   ierr = MatDenseGetArray(St,&array);CHKERRQ(ierr);
181459ac8732SStefano Zampini   if (!mumps->sym) { /* MUMPS always return a full matrix */
18156444a565SStefano Zampini     if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
18166444a565SStefano Zampini       PetscInt i,j,N=mumps->id.size_schur;
18176444a565SStefano Zampini       for (i=0;i<N;i++) {
18186444a565SStefano Zampini         for (j=0;j<N;j++) {
18196444a565SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
18206444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j];
18216444a565SStefano Zampini #else
18226444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
18236444a565SStefano Zampini #endif
18246444a565SStefano Zampini           array[j*N+i] = val;
18256444a565SStefano Zampini         }
18266444a565SStefano Zampini       }
18276444a565SStefano Zampini     } else { /* stored by columns */
18286444a565SStefano Zampini       ierr = PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
18296444a565SStefano Zampini     }
18306444a565SStefano Zampini   } else { /* either full or lower-triangular (not packed) */
18316444a565SStefano Zampini     if (mumps->id.ICNTL(19) == 2) { /* lower triangular stored by columns */
18326444a565SStefano Zampini       PetscInt i,j,N=mumps->id.size_schur;
18336444a565SStefano Zampini       for (i=0;i<N;i++) {
18346444a565SStefano Zampini         for (j=i;j<N;j++) {
18356444a565SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
18366444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j];
18376444a565SStefano Zampini #else
18386444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
18396444a565SStefano Zampini #endif
18406444a565SStefano Zampini           array[i*N+j] = val;
18416444a565SStefano Zampini           array[j*N+i] = val;
18426444a565SStefano Zampini         }
18436444a565SStefano Zampini       }
18446444a565SStefano Zampini     } else if (mumps->id.ICNTL(19) == 3) { /* full matrix */
18456444a565SStefano Zampini       ierr = PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
18466444a565SStefano Zampini     } else { /* ICNTL(19) == 1 lower triangular stored by rows */
18476444a565SStefano Zampini       PetscInt i,j,N=mumps->id.size_schur;
18486444a565SStefano Zampini       for (i=0;i<N;i++) {
18496444a565SStefano Zampini         for (j=0;j<i+1;j++) {
18506444a565SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
18516444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j];
18526444a565SStefano Zampini #else
18536444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
18546444a565SStefano Zampini #endif
18556444a565SStefano Zampini           array[i*N+j] = val;
18566444a565SStefano Zampini           array[j*N+i] = val;
18576444a565SStefano Zampini         }
18586444a565SStefano Zampini       }
18596444a565SStefano Zampini     }
18606444a565SStefano Zampini   }
18616444a565SStefano Zampini   ierr = MatDenseRestoreArray(St,&array);CHKERRQ(ierr);
18626444a565SStefano Zampini   *S = St;
18636444a565SStefano Zampini   PetscFunctionReturn(0);
18646444a565SStefano Zampini }
18656444a565SStefano Zampini 
186659ac8732SStefano Zampini /* -------------------------------------------------------------------------------------------*/
18675a05ddb0SStefano Zampini PetscErrorCode MatFactorGetSchurComplement_MUMPS(Mat F,Mat* S)
186859ac8732SStefano Zampini {
186959ac8732SStefano Zampini   Mat            St;
1870e69c285eSBarry Smith   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
187159ac8732SStefano Zampini   PetscErrorCode ierr;
187259ac8732SStefano Zampini 
187359ac8732SStefano Zampini   PetscFunctionBegin;
18745a05ddb0SStefano Zampini   if (!mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it");
18757bc32933SStefano Zampini   /* It should be the responsibility of the user to handle different ICNTL(19) cases and factorization stages if they want to work with the raw data */
1876241dbb5eSStefano Zampini   ierr = MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.size_schur,(PetscScalar*)mumps->id.schur,&St);CHKERRQ(ierr);
187759ac8732SStefano Zampini   *S = St;
187859ac8732SStefano Zampini   PetscFunctionReturn(0);
187959ac8732SStefano Zampini }
188059ac8732SStefano Zampini 
188159ac8732SStefano Zampini /* -------------------------------------------------------------------------------------------*/
18826dba178dSStefano Zampini PetscErrorCode MatFactorFactorizeSchurComplement_MUMPS(Mat F)
18833280fbe9SStefano Zampini {
18849481b12dSStefano Zampini   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
18853280fbe9SStefano Zampini   PetscErrorCode ierr;
18863280fbe9SStefano Zampini 
18873280fbe9SStefano Zampini   PetscFunctionBegin;
18883280fbe9SStefano Zampini   if (!mumps->id.ICNTL(19)) { /* do nothing */
18893280fbe9SStefano Zampini     PetscFunctionReturn(0);
18903280fbe9SStefano Zampini   }
1891241dbb5eSStefano Zampini   if (mumps->size > 1) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_SUP,"Parallel Schur complement not yet supported from PETSc");
18923280fbe9SStefano Zampini   ierr = MatMumpsFactorSchur_Private(mumps);CHKERRQ(ierr);
18933280fbe9SStefano Zampini   PetscFunctionReturn(0);
18943280fbe9SStefano Zampini }
18953280fbe9SStefano Zampini 
18965a05ddb0SStefano Zampini PetscErrorCode MatFactorInvertSchurComplement_MUMPS(Mat F)
189759ac8732SStefano Zampini {
1898e69c285eSBarry Smith   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
189959ac8732SStefano Zampini   PetscErrorCode ierr;
190059ac8732SStefano Zampini 
190159ac8732SStefano Zampini   PetscFunctionBegin;
190259ac8732SStefano Zampini   if (!mumps->id.ICNTL(19)) { /* do nothing */
190359ac8732SStefano Zampini     PetscFunctionReturn(0);
190459ac8732SStefano Zampini   }
1905241dbb5eSStefano Zampini   if (mumps->size > 1) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_SUP,"Parallel Schur complement not yet supported from PETSc");
190659ac8732SStefano Zampini   ierr = MatMumpsInvertSchur_Private(mumps);CHKERRQ(ierr);
190759ac8732SStefano Zampini   PetscFunctionReturn(0);
190859ac8732SStefano Zampini }
190959ac8732SStefano Zampini 
19106444a565SStefano Zampini /* -------------------------------------------------------------------------------------------*/
19115a05ddb0SStefano Zampini PetscErrorCode MatFactorSolveSchurComplement_MUMPS(Mat F, Vec rhs, Vec sol)
1912e807eca7SStefano Zampini {
1913e69c285eSBarry Smith   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1914e807eca7SStefano Zampini   MumpsScalar    *orhs;
1915e807eca7SStefano Zampini   PetscScalar    *osol,*nrhs,*nsol;
19169a3a5937SStefano Zampini   PetscInt       orhs_size,osol_size,olrhs_size;
1917e807eca7SStefano Zampini   PetscErrorCode ierr;
1918e807eca7SStefano Zampini 
1919e807eca7SStefano Zampini   PetscFunctionBegin;
19205a05ddb0SStefano Zampini   if (!mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it");
1921241dbb5eSStefano Zampini   if (mumps->size > 1) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_SUP,"Parallel Schur complement not yet supported from PETSc");
19226f3cc6f9SBarry Smith 
1923e807eca7SStefano Zampini   /* swap pointers */
1924e807eca7SStefano Zampini   orhs = mumps->id.redrhs;
19259a3a5937SStefano Zampini   olrhs_size = mumps->id.lredrhs;
1926a12f35bfSStefano Zampini   orhs_size = mumps->sizeredrhs;
1927e807eca7SStefano Zampini   osol = mumps->schur_sol;
1928a12f35bfSStefano Zampini   osol_size = mumps->schur_sizesol;
1929e807eca7SStefano Zampini   ierr = VecGetArray(rhs,&nrhs);CHKERRQ(ierr);
1930e807eca7SStefano Zampini   ierr = VecGetArray(sol,&nsol);CHKERRQ(ierr);
1931e807eca7SStefano Zampini   mumps->id.redrhs = (MumpsScalar*)nrhs;
1932a12f35bfSStefano Zampini   ierr = VecGetLocalSize(rhs,&mumps->sizeredrhs);CHKERRQ(ierr);
19339a3a5937SStefano Zampini   mumps->id.lredrhs = mumps->sizeredrhs;
1934e807eca7SStefano Zampini   mumps->schur_sol = nsol;
1935a12f35bfSStefano Zampini   ierr = VecGetLocalSize(sol,&mumps->schur_sizesol);CHKERRQ(ierr);
1936a12f35bfSStefano Zampini 
1937e807eca7SStefano Zampini   /* solve Schur complement */
1938e807eca7SStefano Zampini   mumps->id.nrhs = 1;
1939e807eca7SStefano Zampini   ierr = MatMumpsSolveSchur_Private(mumps,PETSC_FALSE);CHKERRQ(ierr);
1940e807eca7SStefano Zampini   /* restore pointers */
1941e807eca7SStefano Zampini   ierr = VecRestoreArray(rhs,&nrhs);CHKERRQ(ierr);
1942e807eca7SStefano Zampini   ierr = VecRestoreArray(sol,&nsol);CHKERRQ(ierr);
1943e807eca7SStefano Zampini   mumps->id.redrhs = orhs;
19449a3a5937SStefano Zampini   mumps->id.lredrhs = olrhs_size;
1945a12f35bfSStefano Zampini   mumps->sizeredrhs = orhs_size;
1946e807eca7SStefano Zampini   mumps->schur_sol = osol;
1947a12f35bfSStefano Zampini   mumps->schur_sizesol = osol_size;
1948e807eca7SStefano Zampini   PetscFunctionReturn(0);
1949e807eca7SStefano Zampini }
1950e807eca7SStefano Zampini 
19517404bcfbSStefano Zampini /* -------------------------------------------------------------------------------------------*/
19525a05ddb0SStefano Zampini PetscErrorCode MatFactorSolveSchurComplementTranspose_MUMPS(Mat F, Vec rhs, Vec sol)
19537404bcfbSStefano Zampini {
1954e69c285eSBarry Smith   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
19557404bcfbSStefano Zampini   MumpsScalar    *orhs;
19567404bcfbSStefano Zampini   PetscScalar    *osol,*nrhs,*nsol;
1957a12f35bfSStefano Zampini   PetscInt       orhs_size,osol_size;
19587404bcfbSStefano Zampini   PetscErrorCode ierr;
19597404bcfbSStefano Zampini 
19607404bcfbSStefano Zampini   PetscFunctionBegin;
19615a05ddb0SStefano Zampini   if (!mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it");
1962241dbb5eSStefano Zampini   if (mumps->size > 1) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_SUP,"Parallel Schur complement not yet supported from PETSc");
19636f3cc6f9SBarry Smith 
19647404bcfbSStefano Zampini   /* swap pointers */
19657404bcfbSStefano Zampini   orhs = mumps->id.redrhs;
1966a12f35bfSStefano Zampini   orhs_size = mumps->sizeredrhs;
19677404bcfbSStefano Zampini   osol = mumps->schur_sol;
1968a12f35bfSStefano Zampini   osol_size = mumps->schur_sizesol;
19697404bcfbSStefano Zampini   ierr = VecGetArray(rhs,&nrhs);CHKERRQ(ierr);
19707404bcfbSStefano Zampini   ierr = VecGetArray(sol,&nsol);CHKERRQ(ierr);
19717404bcfbSStefano Zampini   mumps->id.redrhs = (MumpsScalar*)nrhs;
1972a12f35bfSStefano Zampini   ierr = VecGetLocalSize(rhs,&mumps->sizeredrhs);CHKERRQ(ierr);
19737404bcfbSStefano Zampini   mumps->schur_sol = nsol;
1974a12f35bfSStefano Zampini   ierr = VecGetLocalSize(sol,&mumps->schur_sizesol);CHKERRQ(ierr);
1975a12f35bfSStefano Zampini 
19767404bcfbSStefano Zampini   /* solve Schur complement */
19777404bcfbSStefano Zampini   mumps->id.nrhs = 1;
19787404bcfbSStefano Zampini   mumps->id.ICNTL(9) = 0;
19797404bcfbSStefano Zampini   ierr = MatMumpsSolveSchur_Private(mumps,PETSC_FALSE);CHKERRQ(ierr);
19807404bcfbSStefano Zampini   mumps->id.ICNTL(9) = 1;
19817404bcfbSStefano Zampini   /* restore pointers */
19827404bcfbSStefano Zampini   ierr = VecRestoreArray(rhs,&nrhs);CHKERRQ(ierr);
19837404bcfbSStefano Zampini   ierr = VecRestoreArray(sol,&nsol);CHKERRQ(ierr);
19847404bcfbSStefano Zampini   mumps->id.redrhs = orhs;
1985a12f35bfSStefano Zampini   mumps->sizeredrhs = orhs_size;
19867404bcfbSStefano Zampini   mumps->schur_sol = osol;
1987a12f35bfSStefano Zampini   mumps->schur_sizesol = osol_size;
19887404bcfbSStefano Zampini   PetscFunctionReturn(0);
19897404bcfbSStefano Zampini }
19907404bcfbSStefano Zampini 
1991a0b0af32SStefano Zampini /* -------------------------------------------------------------------------------------------*/
1992e8ade678SStefano Zampini PetscErrorCode MatFactorSetSchurComplementSolverType_MUMPS(Mat F, PetscInt sym)
1993a0b0af32SStefano Zampini {
19949481b12dSStefano Zampini   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
1995a0b0af32SStefano Zampini 
1996a0b0af32SStefano Zampini   PetscFunctionBegin;
1997a0b0af32SStefano Zampini   if (mumps->schur_factored && mumps->sym != mumps->schur_sym) {
1998e8ade678SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONG,"Cannot change the Schur solver! Schur complement data has been already factored");
1999a0b0af32SStefano Zampini   }
2000a0b0af32SStefano Zampini   mumps->schur_sym = sym;
2001a0b0af32SStefano Zampini   PetscFunctionReturn(0);
2002a0b0af32SStefano Zampini }
2003a0b0af32SStefano Zampini 
2004e807eca7SStefano Zampini /* -------------------------------------------------------------------------------------------*/
20055ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival)
20065ccb76cbSHong Zhang {
2007e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
20085ccb76cbSHong Zhang 
20095ccb76cbSHong Zhang   PetscFunctionBegin;
2010a5e57a09SHong Zhang   mumps->id.ICNTL(icntl) = ival;
20115ccb76cbSHong Zhang   PetscFunctionReturn(0);
20125ccb76cbSHong Zhang }
20135ccb76cbSHong Zhang 
2014bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt *ival)
2015bc6112feSHong Zhang {
2016e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
2017bc6112feSHong Zhang 
2018bc6112feSHong Zhang   PetscFunctionBegin;
2019bc6112feSHong Zhang   *ival = mumps->id.ICNTL(icntl);
2020bc6112feSHong Zhang   PetscFunctionReturn(0);
2021bc6112feSHong Zhang }
2022bc6112feSHong Zhang 
20235ccb76cbSHong Zhang /*@
20245ccb76cbSHong Zhang   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()
20255ccb76cbSHong Zhang 
20265ccb76cbSHong Zhang    Logically Collective on Mat
20275ccb76cbSHong Zhang 
20285ccb76cbSHong Zhang    Input Parameters:
20295ccb76cbSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
20305ccb76cbSHong Zhang .  icntl - index of MUMPS parameter array ICNTL()
20315ccb76cbSHong Zhang -  ival - value of MUMPS ICNTL(icntl)
20325ccb76cbSHong Zhang 
20335ccb76cbSHong Zhang   Options Database:
20345ccb76cbSHong Zhang .   -mat_mumps_icntl_<icntl> <ival>
20355ccb76cbSHong Zhang 
20365ccb76cbSHong Zhang    Level: beginner
20375ccb76cbSHong Zhang 
203896a0c994SBarry Smith    References:
203996a0c994SBarry Smith .     MUMPS Users' Guide
20405ccb76cbSHong Zhang 
20419fc87aa7SBarry Smith .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
20425ccb76cbSHong Zhang  @*/
20435ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival)
20445ccb76cbSHong Zhang {
20455ccb76cbSHong Zhang   PetscErrorCode ierr;
20465ccb76cbSHong Zhang 
20475ccb76cbSHong Zhang   PetscFunctionBegin;
20482989dfd4SHong Zhang   PetscValidType(F,1);
20492989dfd4SHong Zhang   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
20505ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
20515ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,ival,3);
20525ccb76cbSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr);
20535ccb76cbSHong Zhang   PetscFunctionReturn(0);
20545ccb76cbSHong Zhang }
20555ccb76cbSHong Zhang 
2056a21f80fcSHong Zhang /*@
2057a21f80fcSHong Zhang   MatMumpsGetIcntl - Get MUMPS parameter ICNTL()
2058a21f80fcSHong Zhang 
2059a21f80fcSHong Zhang    Logically Collective on Mat
2060a21f80fcSHong Zhang 
2061a21f80fcSHong Zhang    Input Parameters:
2062a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2063a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array ICNTL()
2064a21f80fcSHong Zhang 
2065a21f80fcSHong Zhang   Output Parameter:
2066a21f80fcSHong Zhang .  ival - value of MUMPS ICNTL(icntl)
2067a21f80fcSHong Zhang 
2068a21f80fcSHong Zhang    Level: beginner
2069a21f80fcSHong Zhang 
207096a0c994SBarry Smith    References:
207196a0c994SBarry Smith .     MUMPS Users' Guide
2072a21f80fcSHong Zhang 
20739fc87aa7SBarry Smith .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2074a21f80fcSHong Zhang @*/
2075bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl(Mat F,PetscInt icntl,PetscInt *ival)
2076bc6112feSHong Zhang {
2077bc6112feSHong Zhang   PetscErrorCode ierr;
2078bc6112feSHong Zhang 
2079bc6112feSHong Zhang   PetscFunctionBegin;
20802989dfd4SHong Zhang   PetscValidType(F,1);
20812989dfd4SHong Zhang   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2082bc6112feSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
2083bc6112feSHong Zhang   PetscValidIntPointer(ival,3);
20842989dfd4SHong Zhang   ierr = PetscUseMethod(F,"MatMumpsGetIcntl_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
2085bc6112feSHong Zhang   PetscFunctionReturn(0);
2086bc6112feSHong Zhang }
2087bc6112feSHong Zhang 
20888928b65cSHong Zhang /* -------------------------------------------------------------------------------------------*/
20898928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal val)
20908928b65cSHong Zhang {
2091e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
20928928b65cSHong Zhang 
20938928b65cSHong Zhang   PetscFunctionBegin;
20948928b65cSHong Zhang   mumps->id.CNTL(icntl) = val;
20958928b65cSHong Zhang   PetscFunctionReturn(0);
20968928b65cSHong Zhang }
20978928b65cSHong Zhang 
2098bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal *val)
2099bc6112feSHong Zhang {
2100e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
2101bc6112feSHong Zhang 
2102bc6112feSHong Zhang   PetscFunctionBegin;
2103bc6112feSHong Zhang   *val = mumps->id.CNTL(icntl);
2104bc6112feSHong Zhang   PetscFunctionReturn(0);
2105bc6112feSHong Zhang }
2106bc6112feSHong Zhang 
21078928b65cSHong Zhang /*@
21088928b65cSHong Zhang   MatMumpsSetCntl - Set MUMPS parameter CNTL()
21098928b65cSHong Zhang 
21108928b65cSHong Zhang    Logically Collective on Mat
21118928b65cSHong Zhang 
21128928b65cSHong Zhang    Input Parameters:
21138928b65cSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
21148928b65cSHong Zhang .  icntl - index of MUMPS parameter array CNTL()
21158928b65cSHong Zhang -  val - value of MUMPS CNTL(icntl)
21168928b65cSHong Zhang 
21178928b65cSHong Zhang   Options Database:
21188928b65cSHong Zhang .   -mat_mumps_cntl_<icntl> <val>
21198928b65cSHong Zhang 
21208928b65cSHong Zhang    Level: beginner
21218928b65cSHong Zhang 
212296a0c994SBarry Smith    References:
212396a0c994SBarry Smith .     MUMPS Users' Guide
21248928b65cSHong Zhang 
21259fc87aa7SBarry Smith .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
21268928b65cSHong Zhang @*/
21278928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl(Mat F,PetscInt icntl,PetscReal val)
21288928b65cSHong Zhang {
21298928b65cSHong Zhang   PetscErrorCode ierr;
21308928b65cSHong Zhang 
21318928b65cSHong Zhang   PetscFunctionBegin;
21322989dfd4SHong Zhang   PetscValidType(F,1);
21332989dfd4SHong Zhang   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
21348928b65cSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
2135bc6112feSHong Zhang   PetscValidLogicalCollectiveReal(F,val,3);
21368928b65cSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsSetCntl_C",(Mat,PetscInt,PetscReal),(F,icntl,val));CHKERRQ(ierr);
21378928b65cSHong Zhang   PetscFunctionReturn(0);
21388928b65cSHong Zhang }
21398928b65cSHong Zhang 
2140a21f80fcSHong Zhang /*@
2141a21f80fcSHong Zhang   MatMumpsGetCntl - Get MUMPS parameter CNTL()
2142a21f80fcSHong Zhang 
2143a21f80fcSHong Zhang    Logically Collective on Mat
2144a21f80fcSHong Zhang 
2145a21f80fcSHong Zhang    Input Parameters:
2146a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2147a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array CNTL()
2148a21f80fcSHong Zhang 
2149a21f80fcSHong Zhang   Output Parameter:
2150a21f80fcSHong Zhang .  val - value of MUMPS CNTL(icntl)
2151a21f80fcSHong Zhang 
2152a21f80fcSHong Zhang    Level: beginner
2153a21f80fcSHong Zhang 
215496a0c994SBarry Smith    References:
215596a0c994SBarry Smith .      MUMPS Users' Guide
2156a21f80fcSHong Zhang 
21579fc87aa7SBarry Smith .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2158a21f80fcSHong Zhang @*/
2159bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl(Mat F,PetscInt icntl,PetscReal *val)
2160bc6112feSHong Zhang {
2161bc6112feSHong Zhang   PetscErrorCode ierr;
2162bc6112feSHong Zhang 
2163bc6112feSHong Zhang   PetscFunctionBegin;
21642989dfd4SHong Zhang   PetscValidType(F,1);
21652989dfd4SHong Zhang   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2166bc6112feSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
2167bc6112feSHong Zhang   PetscValidRealPointer(val,3);
21682989dfd4SHong Zhang   ierr = PetscUseMethod(F,"MatMumpsGetCntl_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
2169bc6112feSHong Zhang   PetscFunctionReturn(0);
2170bc6112feSHong Zhang }
2171bc6112feSHong Zhang 
2172ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F,PetscInt icntl,PetscInt *info)
2173bc6112feSHong Zhang {
2174e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
2175bc6112feSHong Zhang 
2176bc6112feSHong Zhang   PetscFunctionBegin;
2177bc6112feSHong Zhang   *info = mumps->id.INFO(icntl);
2178bc6112feSHong Zhang   PetscFunctionReturn(0);
2179bc6112feSHong Zhang }
2180bc6112feSHong Zhang 
2181ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F,PetscInt icntl,PetscInt *infog)
2182bc6112feSHong Zhang {
2183e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
2184bc6112feSHong Zhang 
2185bc6112feSHong Zhang   PetscFunctionBegin;
2186bc6112feSHong Zhang   *infog = mumps->id.INFOG(icntl);
2187bc6112feSHong Zhang   PetscFunctionReturn(0);
2188bc6112feSHong Zhang }
2189bc6112feSHong Zhang 
2190ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfo)
2191bc6112feSHong Zhang {
2192e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
2193bc6112feSHong Zhang 
2194bc6112feSHong Zhang   PetscFunctionBegin;
2195bc6112feSHong Zhang   *rinfo = mumps->id.RINFO(icntl);
2196bc6112feSHong Zhang   PetscFunctionReturn(0);
2197bc6112feSHong Zhang }
2198bc6112feSHong Zhang 
2199ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfog)
2200bc6112feSHong Zhang {
2201e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
2202bc6112feSHong Zhang 
2203bc6112feSHong Zhang   PetscFunctionBegin;
2204bc6112feSHong Zhang   *rinfog = mumps->id.RINFOG(icntl);
2205bc6112feSHong Zhang   PetscFunctionReturn(0);
2206bc6112feSHong Zhang }
2207bc6112feSHong Zhang 
2208a21f80fcSHong Zhang /*@
2209a21f80fcSHong Zhang   MatMumpsGetInfo - Get MUMPS parameter INFO()
2210a21f80fcSHong Zhang 
2211a21f80fcSHong Zhang    Logically Collective on Mat
2212a21f80fcSHong Zhang 
2213a21f80fcSHong Zhang    Input Parameters:
2214a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2215a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array INFO()
2216a21f80fcSHong Zhang 
2217a21f80fcSHong Zhang   Output Parameter:
2218a21f80fcSHong Zhang .  ival - value of MUMPS INFO(icntl)
2219a21f80fcSHong Zhang 
2220a21f80fcSHong Zhang    Level: beginner
2221a21f80fcSHong Zhang 
222296a0c994SBarry Smith    References:
222396a0c994SBarry Smith .      MUMPS Users' Guide
2224a21f80fcSHong Zhang 
22259fc87aa7SBarry Smith .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2226a21f80fcSHong Zhang @*/
2227ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo(Mat F,PetscInt icntl,PetscInt *ival)
2228bc6112feSHong Zhang {
2229bc6112feSHong Zhang   PetscErrorCode ierr;
2230bc6112feSHong Zhang 
2231bc6112feSHong Zhang   PetscFunctionBegin;
22322989dfd4SHong Zhang   PetscValidType(F,1);
22332989dfd4SHong Zhang   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2234ca810319SHong Zhang   PetscValidIntPointer(ival,3);
22352989dfd4SHong Zhang   ierr = PetscUseMethod(F,"MatMumpsGetInfo_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
2236bc6112feSHong Zhang   PetscFunctionReturn(0);
2237bc6112feSHong Zhang }
2238bc6112feSHong Zhang 
2239a21f80fcSHong Zhang /*@
2240a21f80fcSHong Zhang   MatMumpsGetInfog - Get MUMPS parameter INFOG()
2241a21f80fcSHong Zhang 
2242a21f80fcSHong Zhang    Logically Collective on Mat
2243a21f80fcSHong Zhang 
2244a21f80fcSHong Zhang    Input Parameters:
2245a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2246a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array INFOG()
2247a21f80fcSHong Zhang 
2248a21f80fcSHong Zhang   Output Parameter:
2249a21f80fcSHong Zhang .  ival - value of MUMPS INFOG(icntl)
2250a21f80fcSHong Zhang 
2251a21f80fcSHong Zhang    Level: beginner
2252a21f80fcSHong Zhang 
225396a0c994SBarry Smith    References:
225496a0c994SBarry Smith .      MUMPS Users' Guide
2255a21f80fcSHong Zhang 
22569fc87aa7SBarry Smith .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2257a21f80fcSHong Zhang @*/
2258ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog(Mat F,PetscInt icntl,PetscInt *ival)
2259bc6112feSHong Zhang {
2260bc6112feSHong Zhang   PetscErrorCode ierr;
2261bc6112feSHong Zhang 
2262bc6112feSHong Zhang   PetscFunctionBegin;
22632989dfd4SHong Zhang   PetscValidType(F,1);
22642989dfd4SHong Zhang   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2265ca810319SHong Zhang   PetscValidIntPointer(ival,3);
22662989dfd4SHong Zhang   ierr = PetscUseMethod(F,"MatMumpsGetInfog_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
2267bc6112feSHong Zhang   PetscFunctionReturn(0);
2268bc6112feSHong Zhang }
2269bc6112feSHong Zhang 
2270a21f80fcSHong Zhang /*@
2271a21f80fcSHong Zhang   MatMumpsGetRinfo - Get MUMPS parameter RINFO()
2272a21f80fcSHong Zhang 
2273a21f80fcSHong Zhang    Logically Collective on Mat
2274a21f80fcSHong Zhang 
2275a21f80fcSHong Zhang    Input Parameters:
2276a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2277a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array RINFO()
2278a21f80fcSHong Zhang 
2279a21f80fcSHong Zhang   Output Parameter:
2280a21f80fcSHong Zhang .  val - value of MUMPS RINFO(icntl)
2281a21f80fcSHong Zhang 
2282a21f80fcSHong Zhang    Level: beginner
2283a21f80fcSHong Zhang 
228496a0c994SBarry Smith    References:
228596a0c994SBarry Smith .       MUMPS Users' Guide
2286a21f80fcSHong Zhang 
22879fc87aa7SBarry Smith .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2288a21f80fcSHong Zhang @*/
2289ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo(Mat F,PetscInt icntl,PetscReal *val)
2290bc6112feSHong Zhang {
2291bc6112feSHong Zhang   PetscErrorCode ierr;
2292bc6112feSHong Zhang 
2293bc6112feSHong Zhang   PetscFunctionBegin;
22942989dfd4SHong Zhang   PetscValidType(F,1);
22952989dfd4SHong Zhang   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2296bc6112feSHong Zhang   PetscValidRealPointer(val,3);
22972989dfd4SHong Zhang   ierr = PetscUseMethod(F,"MatMumpsGetRinfo_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
2298bc6112feSHong Zhang   PetscFunctionReturn(0);
2299bc6112feSHong Zhang }
2300bc6112feSHong Zhang 
2301a21f80fcSHong Zhang /*@
2302a21f80fcSHong Zhang   MatMumpsGetRinfog - Get MUMPS parameter RINFOG()
2303a21f80fcSHong Zhang 
2304a21f80fcSHong Zhang    Logically Collective on Mat
2305a21f80fcSHong Zhang 
2306a21f80fcSHong Zhang    Input Parameters:
2307a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2308a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array RINFOG()
2309a21f80fcSHong Zhang 
2310a21f80fcSHong Zhang   Output Parameter:
2311a21f80fcSHong Zhang .  val - value of MUMPS RINFOG(icntl)
2312a21f80fcSHong Zhang 
2313a21f80fcSHong Zhang    Level: beginner
2314a21f80fcSHong Zhang 
231596a0c994SBarry Smith    References:
231696a0c994SBarry Smith .      MUMPS Users' Guide
2317a21f80fcSHong Zhang 
23189fc87aa7SBarry Smith .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2319a21f80fcSHong Zhang @*/
2320ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog(Mat F,PetscInt icntl,PetscReal *val)
2321bc6112feSHong Zhang {
2322bc6112feSHong Zhang   PetscErrorCode ierr;
2323bc6112feSHong Zhang 
2324bc6112feSHong Zhang   PetscFunctionBegin;
23252989dfd4SHong Zhang   PetscValidType(F,1);
23262989dfd4SHong Zhang   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2327bc6112feSHong Zhang   PetscValidRealPointer(val,3);
23282989dfd4SHong Zhang   ierr = PetscUseMethod(F,"MatMumpsGetRinfog_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
2329bc6112feSHong Zhang   PetscFunctionReturn(0);
2330bc6112feSHong Zhang }
2331bc6112feSHong Zhang 
233224b6179bSKris Buschelman /*MC
23332692d6eeSBarry Smith   MATSOLVERMUMPS -  A matrix type providing direct solvers (LU and Cholesky) for
233424b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
233524b6179bSKris Buschelman 
233641c8de11SBarry Smith   Works with MATAIJ and MATSBAIJ matrices
233724b6179bSKris Buschelman 
2338c2b89b5dSBarry Smith   Use ./configure --download-mumps --download-scalapack --download-parmetis --download-metis --download-ptscotch  to have PETSc installed with MUMPS
2339c2b89b5dSBarry Smith 
2340c2b89b5dSBarry Smith   Use -pc_type cholesky or lu -pc_factor_mat_solver_package mumps to us this direct solver
2341c2b89b5dSBarry Smith 
234224b6179bSKris Buschelman   Options Database Keys:
23434422a9fcSPatrick Sanan +  -mat_mumps_icntl_1 - ICNTL(1): output stream for error messages
23444422a9fcSPatrick Sanan .  -mat_mumps_icntl_2 - ICNTL(2): output stream for diagnostic printing, statistics, and warning
23454422a9fcSPatrick Sanan .  -mat_mumps_icntl_3 -  ICNTL(3): output stream for global information, collected on the host
23464422a9fcSPatrick Sanan .  -mat_mumps_icntl_4 -  ICNTL(4): level of printing (0 to 4)
23474422a9fcSPatrick Sanan .  -mat_mumps_icntl_6 - ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7)
23484422a9fcSPatrick Sanan .  -mat_mumps_icntl_7 - ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis
23494422a9fcSPatrick Sanan .  -mat_mumps_icntl_8  - ICNTL(8): scaling strategy (-2 to 8 or 77)
23504422a9fcSPatrick Sanan .  -mat_mumps_icntl_10  - ICNTL(10): max num of refinements
23514422a9fcSPatrick Sanan .  -mat_mumps_icntl_11  - ICNTL(11): statistics related to an error analysis (via -ksp_view)
23524422a9fcSPatrick Sanan .  -mat_mumps_icntl_12  - ICNTL(12): an ordering strategy for symmetric matrices (0 to 3)
23534422a9fcSPatrick Sanan .  -mat_mumps_icntl_13  - ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting
23544422a9fcSPatrick Sanan .  -mat_mumps_icntl_14  - ICNTL(14): percentage increase in the estimated working space
23554422a9fcSPatrick Sanan .  -mat_mumps_icntl_19  - ICNTL(19): computes the Schur complement
23564422a9fcSPatrick Sanan .  -mat_mumps_icntl_22  - ICNTL(22): in-core/out-of-core factorization and solve (0 or 1)
23574422a9fcSPatrick Sanan .  -mat_mumps_icntl_23  - ICNTL(23): max size of the working memory (MB) that can allocate per processor
23584422a9fcSPatrick Sanan .  -mat_mumps_icntl_24  - ICNTL(24): detection of null pivot rows (0 or 1)
23594422a9fcSPatrick Sanan .  -mat_mumps_icntl_25  - ICNTL(25): compute a solution of a deficient matrix and a null space basis
23604422a9fcSPatrick Sanan .  -mat_mumps_icntl_26  - ICNTL(26): drives the solution phase if a Schur complement matrix
23614422a9fcSPatrick Sanan .  -mat_mumps_icntl_28  - ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering
23624422a9fcSPatrick Sanan .  -mat_mumps_icntl_29 - ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis
23634422a9fcSPatrick Sanan .  -mat_mumps_icntl_30 - ICNTL(30): compute user-specified set of entries in inv(A)
23644422a9fcSPatrick Sanan .  -mat_mumps_icntl_31 - ICNTL(31): indicates which factors may be discarded during factorization
23654422a9fcSPatrick Sanan .  -mat_mumps_icntl_33 - ICNTL(33): compute determinant
23664422a9fcSPatrick Sanan .  -mat_mumps_cntl_1  - CNTL(1): relative pivoting threshold
23674422a9fcSPatrick Sanan .  -mat_mumps_cntl_2  -  CNTL(2): stopping criterion of refinement
23684422a9fcSPatrick Sanan .  -mat_mumps_cntl_3 - CNTL(3): absolute pivoting threshold
23694422a9fcSPatrick Sanan .  -mat_mumps_cntl_4 - CNTL(4): value for static pivoting
23704422a9fcSPatrick Sanan -  -mat_mumps_cntl_5 - CNTL(5): fixation for null pivots
237124b6179bSKris Buschelman 
237224b6179bSKris Buschelman   Level: beginner
237324b6179bSKris Buschelman 
23749fc87aa7SBarry Smith     Notes: When a MUMPS factorization fails inside a KSP solve, for example with a KSP_DIVERGED_PCSETUP_FAILED, one can find the MUMPS information about the failure by calling
23759fc87aa7SBarry Smith $          KSPGetPC(ksp,&pc);
23769fc87aa7SBarry Smith $          PCFactorGetMatrix(pc,&mat);
23779fc87aa7SBarry Smith $          MatMumpsGetInfo(mat,....);
23789fc87aa7SBarry Smith $          MatMumpsGetInfog(mat,....); etc.
23799fc87aa7SBarry Smith            Or you can run with -ksp_error_if_not_converged and the program will be stopped and the information printed in the error message.
23809fc87aa7SBarry Smith 
23819fc87aa7SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage, MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog(), KSPGetPC(), PCGetFactor(), PCFactorGetMatrix()
238241c8de11SBarry Smith 
238324b6179bSKris Buschelman M*/
238424b6179bSKris Buschelman 
2385f7a08781SBarry Smith static PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type)
238635bd34faSBarry Smith {
238735bd34faSBarry Smith   PetscFunctionBegin;
23882692d6eeSBarry Smith   *type = MATSOLVERMUMPS;
238935bd34faSBarry Smith   PetscFunctionReturn(0);
239035bd34faSBarry Smith }
239135bd34faSBarry Smith 
2392bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI AIJ matrices */
2393cc2e6a90SBarry Smith static PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
23942877fffaSHong Zhang {
23952877fffaSHong Zhang   Mat            B;
23962877fffaSHong Zhang   PetscErrorCode ierr;
23972877fffaSHong Zhang   Mat_MUMPS      *mumps;
2398ace3abfcSBarry Smith   PetscBool      isSeqAIJ;
23992877fffaSHong Zhang 
24002877fffaSHong Zhang   PetscFunctionBegin;
24012877fffaSHong Zhang   /* Create the factorization matrix */
2402251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr);
2403ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
24042877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
2405e69c285eSBarry Smith   ierr = PetscStrallocpy("mumps",&((PetscObject)B)->type_name);CHKERRQ(ierr);
2406e69c285eSBarry Smith   ierr = MatSetUp(B);CHKERRQ(ierr);
24072877fffaSHong Zhang 
2408b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
24092205254eSKarl Rupp 
24102877fffaSHong Zhang   B->ops->view        = MatView_MUMPS;
241135bd34faSBarry Smith   B->ops->getinfo     = MatGetInfo_MUMPS;
24122205254eSKarl Rupp 
2413bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
24145a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr);
24155a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorInvertSchurComplement_C",MatFactorInvertSchurComplement_MUMPS);CHKERRQ(ierr);
24165a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr);
24175a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSchurComplement_C",MatFactorGetSchurComplement_MUMPS);CHKERRQ(ierr);
24185a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplement_C",MatFactorSolveSchurComplement_MUMPS);CHKERRQ(ierr);
24195a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplementTranspose_C",MatFactorSolveSchurComplementTranspose_MUMPS);CHKERRQ(ierr);
24206dba178dSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorFactorizeSchurComplement_C",MatFactorFactorizeSchurComplement_MUMPS);CHKERRQ(ierr);
2421e8ade678SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurComplementSolverType_C",MatFactorSetSchurComplementSolverType_MUMPS);CHKERRQ(ierr);
2422bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
2423bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
2424bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
2425bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
2426ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
2427ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
2428ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
2429ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
24306444a565SStefano Zampini 
2431450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
2432450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
2433d5f3da31SBarry Smith     B->factortype            = MAT_FACTOR_LU;
2434bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
2435bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
2436746480a1SHong Zhang     mumps->sym = 0;
2437dcd589f8SShri Abhyankar   } else {
243867877ebaSShri Abhyankar     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2439450b117fSShri Abhyankar     B->factortype                  = MAT_FACTOR_CHOLESKY;
2440bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
2441bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
244259ac8732SStefano Zampini #if defined(PETSC_USE_COMPLEX)
244359ac8732SStefano Zampini     mumps->sym = 2;
244459ac8732SStefano Zampini #else
24456fdc2a6dSBarry Smith     if (A->spd_set && A->spd) mumps->sym = 1;
24466fdc2a6dSBarry Smith     else                      mumps->sym = 2;
244759ac8732SStefano Zampini #endif
2448450b117fSShri Abhyankar   }
24492877fffaSHong Zhang 
245000c67f3bSHong Zhang   /* set solvertype */
245100c67f3bSHong Zhang   ierr = PetscFree(B->solvertype);CHKERRQ(ierr);
245200c67f3bSHong Zhang   ierr = PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);CHKERRQ(ierr);
245300c67f3bSHong Zhang 
24542877fffaSHong Zhang   mumps->isAIJ    = PETSC_TRUE;
24552877fffaSHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
2456e69c285eSBarry Smith   B->data        = (void*)mumps;
24572205254eSKarl Rupp 
2458f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
2459746480a1SHong Zhang 
24602877fffaSHong Zhang   *F = B;
24612877fffaSHong Zhang   PetscFunctionReturn(0);
24622877fffaSHong Zhang }
24632877fffaSHong Zhang 
2464bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI SBAIJ matrices */
2465cc2e6a90SBarry Smith static PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
24662877fffaSHong Zhang {
24672877fffaSHong Zhang   Mat            B;
24682877fffaSHong Zhang   PetscErrorCode ierr;
24692877fffaSHong Zhang   Mat_MUMPS      *mumps;
2470ace3abfcSBarry Smith   PetscBool      isSeqSBAIJ;
24712877fffaSHong Zhang 
24722877fffaSHong Zhang   PetscFunctionBegin;
2473ce94432eSBarry Smith   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
2474ce94432eSBarry Smith   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");
2475251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
24762877fffaSHong Zhang   /* Create the factorization matrix */
2477ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
24782877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
2479e69c285eSBarry Smith   ierr = PetscStrallocpy("mumps",&((PetscObject)B)->type_name);CHKERRQ(ierr);
2480e69c285eSBarry Smith   ierr = MatSetUp(B);CHKERRQ(ierr);
2481e69c285eSBarry Smith 
2482b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
2483bccb9932SShri Abhyankar   if (isSeqSBAIJ) {
248416ebf90aSShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
2485dcd589f8SShri Abhyankar   } else {
2486bccb9932SShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
2487bccb9932SShri Abhyankar   }
2488bccb9932SShri Abhyankar 
2489e69c285eSBarry Smith   B->ops->getinfo                = MatGetInfo_External;
249067877ebaSShri Abhyankar   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2491bccb9932SShri Abhyankar   B->ops->view                   = MatView_MUMPS;
24922205254eSKarl Rupp 
2493bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
24945a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr);
24955a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorInvertSchurComplement_C",MatFactorInvertSchurComplement_MUMPS);CHKERRQ(ierr);
24965a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr);
24975a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSchurComplement_C",MatFactorGetSchurComplement_MUMPS);CHKERRQ(ierr);
24985a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplement_C",MatFactorSolveSchurComplement_MUMPS);CHKERRQ(ierr);
24995a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplementTranspose_C",MatFactorSolveSchurComplementTranspose_MUMPS);CHKERRQ(ierr);
25006dba178dSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorFactorizeSchurComplement_C",MatFactorFactorizeSchurComplement_MUMPS);CHKERRQ(ierr);
2501e8ade678SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurComplementSolverType_C",MatFactorSetSchurComplementSolverType_MUMPS);CHKERRQ(ierr);
2502b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
2503b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
2504b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
2505b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
2506ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
2507ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
2508ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
2509ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
25102205254eSKarl Rupp 
2511f4762488SHong Zhang   B->factortype = MAT_FACTOR_CHOLESKY;
251259ac8732SStefano Zampini #if defined(PETSC_USE_COMPLEX)
251359ac8732SStefano Zampini   mumps->sym = 2;
251459ac8732SStefano Zampini #else
25156fdc2a6dSBarry Smith   if (A->spd_set && A->spd) mumps->sym = 1;
25166fdc2a6dSBarry Smith   else                      mumps->sym = 2;
251759ac8732SStefano Zampini #endif
2518a214ac2aSShri Abhyankar 
251900c67f3bSHong Zhang   /* set solvertype */
252000c67f3bSHong Zhang   ierr = PetscFree(B->solvertype);CHKERRQ(ierr);
252100c67f3bSHong Zhang   ierr = PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);CHKERRQ(ierr);
252200c67f3bSHong Zhang 
2523bccb9932SShri Abhyankar   mumps->isAIJ    = PETSC_FALSE;
2524f3c0ef26SHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
2525e69c285eSBarry Smith   B->data        = (void*)mumps;
25262205254eSKarl Rupp 
2527f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
2528746480a1SHong Zhang 
25292877fffaSHong Zhang   *F = B;
25302877fffaSHong Zhang   PetscFunctionReturn(0);
25312877fffaSHong Zhang }
253297969023SHong Zhang 
2533cc2e6a90SBarry Smith static PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
253467877ebaSShri Abhyankar {
253567877ebaSShri Abhyankar   Mat            B;
253667877ebaSShri Abhyankar   PetscErrorCode ierr;
253767877ebaSShri Abhyankar   Mat_MUMPS      *mumps;
2538ace3abfcSBarry Smith   PetscBool      isSeqBAIJ;
253967877ebaSShri Abhyankar 
254067877ebaSShri Abhyankar   PetscFunctionBegin;
254167877ebaSShri Abhyankar   /* Create the factorization matrix */
2542251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr);
2543ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
254467877ebaSShri Abhyankar   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
2545e69c285eSBarry Smith   ierr = PetscStrallocpy("mumps",&((PetscObject)B)->type_name);CHKERRQ(ierr);
2546e69c285eSBarry Smith   ierr = MatSetUp(B);CHKERRQ(ierr);
2547450b117fSShri Abhyankar 
2548b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
2549450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
2550450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
2551450b117fSShri Abhyankar     B->factortype            = MAT_FACTOR_LU;
2552bccb9932SShri Abhyankar     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
2553bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
2554746480a1SHong Zhang     mumps->sym = 0;
2555f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n");
2556bccb9932SShri Abhyankar 
2557e69c285eSBarry Smith   B->ops->getinfo     = MatGetInfo_External;
2558450b117fSShri Abhyankar   B->ops->view        = MatView_MUMPS;
25592205254eSKarl Rupp 
2560bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
25615a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr);
25625a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorInvertSchurComplement_C",MatFactorInvertSchurComplement_MUMPS);CHKERRQ(ierr);
25635a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr);
25645a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSchurComplement_C",MatFactorGetSchurComplement_MUMPS);CHKERRQ(ierr);
25655a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplement_C",MatFactorSolveSchurComplement_MUMPS);CHKERRQ(ierr);
25665a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplementTranspose_C",MatFactorSolveSchurComplementTranspose_MUMPS);CHKERRQ(ierr);
25676dba178dSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorFactorizeSchurComplement_C",MatFactorFactorizeSchurComplement_MUMPS);CHKERRQ(ierr);
2568e8ade678SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurComplementSolverType_C",MatFactorSetSchurComplementSolverType_MUMPS);CHKERRQ(ierr);
2569bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
2570bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
2571bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
2572bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
2573ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
2574ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
2575ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
2576ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
2577450b117fSShri Abhyankar 
257800c67f3bSHong Zhang   /* set solvertype */
257900c67f3bSHong Zhang   ierr = PetscFree(B->solvertype);CHKERRQ(ierr);
258000c67f3bSHong Zhang   ierr = PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);CHKERRQ(ierr);
258100c67f3bSHong Zhang 
2582450b117fSShri Abhyankar   mumps->isAIJ    = PETSC_TRUE;
2583450b117fSShri Abhyankar   B->ops->destroy = MatDestroy_MUMPS;
2584e69c285eSBarry Smith   B->data        = (void*)mumps;
25852205254eSKarl Rupp 
2586f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
2587746480a1SHong Zhang 
2588450b117fSShri Abhyankar   *F = B;
2589450b117fSShri Abhyankar   PetscFunctionReturn(0);
2590450b117fSShri Abhyankar }
259142c9c57cSBarry Smith 
259229b38603SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_MUMPS(void)
259342c9c57cSBarry Smith {
259442c9c57cSBarry Smith   PetscErrorCode ierr;
259542c9c57cSBarry Smith 
259642c9c57cSBarry Smith   PetscFunctionBegin;
259742c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr);
259842c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr);
259942c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr);
260042c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr);
260142c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPISBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr);
260242c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr);
260342c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr);
260442c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr);
260542c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr);
260642c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr);
260742c9c57cSBarry Smith   PetscFunctionReturn(0);
260842c9c57cSBarry Smith }
260942c9c57cSBarry Smith 
2610