xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision 6dba178d53d088f8af27c39b86adaa0496f0d2e2)
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 
100bf0cc555SLisandro Dalcin   PetscErrorCode (*Destroy)(Mat);
101bccb9932SShri Abhyankar   PetscErrorCode (*ConvertToTriples)(Mat, int, MatReuse, int*, int**, int**, PetscScalar**);
102f0c56d0fSKris Buschelman } Mat_MUMPS;
103f0c56d0fSKris Buschelman 
10409573ac7SBarry Smith extern PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);
105b24902e0SBarry Smith 
10659ac8732SStefano Zampini #undef __FUNCT__
10759ac8732SStefano Zampini #define __FUNCT__ "MatMumpsResetSchur_Private"
10859ac8732SStefano Zampini static PetscErrorCode MatMumpsResetSchur_Private(Mat_MUMPS* mumps)
109b5fa320bSStefano Zampini {
110b5fa320bSStefano Zampini   PetscErrorCode ierr;
111b5fa320bSStefano Zampini 
112b5fa320bSStefano Zampini   PetscFunctionBegin;
11359ac8732SStefano Zampini   ierr = PetscFree2(mumps->id.listvar_schur,mumps->id.schur);CHKERRQ(ierr);
11459ac8732SStefano Zampini   ierr = PetscFree(mumps->id.redrhs);CHKERRQ(ierr);
11559ac8732SStefano Zampini   ierr = PetscFree(mumps->schur_sol);CHKERRQ(ierr);
11659ac8732SStefano Zampini   ierr = PetscFree(mumps->schur_pivots);CHKERRQ(ierr);
11759ac8732SStefano Zampini   ierr = PetscFree(mumps->schur_work);CHKERRQ(ierr);
11859ac8732SStefano Zampini   mumps->id.size_schur = 0;
11959ac8732SStefano Zampini   mumps->id.ICNTL(19) = 0;
12059ac8732SStefano Zampini   PetscFunctionReturn(0);
12159ac8732SStefano Zampini }
12259ac8732SStefano Zampini 
12359ac8732SStefano Zampini #undef __FUNCT__
12459ac8732SStefano Zampini #define __FUNCT__ "MatMumpsFactorSchur_Private"
12559ac8732SStefano Zampini static PetscErrorCode MatMumpsFactorSchur_Private(Mat_MUMPS* mumps)
12659ac8732SStefano Zampini {
12759ac8732SStefano Zampini   PetscBLASInt   B_N,B_ierr,B_slda;
12859ac8732SStefano Zampini   PetscErrorCode ierr;
12959ac8732SStefano Zampini 
13059ac8732SStefano Zampini   PetscFunctionBegin;
13159ac8732SStefano Zampini   if (mumps->schur_factored) {
13259ac8732SStefano Zampini     PetscFunctionReturn(0);
13359ac8732SStefano Zampini   }
13459ac8732SStefano Zampini   ierr = PetscBLASIntCast(mumps->id.size_schur,&B_N);CHKERRQ(ierr);
13559ac8732SStefano Zampini   ierr = PetscBLASIntCast(mumps->id.schur_lld,&B_slda);CHKERRQ(ierr);
13659ac8732SStefano Zampini   if (!mumps->sym) { /* MUMPS always return a full Schur matrix */
13759ac8732SStefano Zampini     if (!mumps->schur_pivots) {
13859ac8732SStefano Zampini       ierr = PetscMalloc1(B_N,&mumps->schur_pivots);CHKERRQ(ierr);
13959ac8732SStefano Zampini     }
14059ac8732SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
14159ac8732SStefano Zampini     PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,(PetscScalar*)mumps->id.schur,&B_slda,mumps->schur_pivots,&B_ierr));
14259ac8732SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
14359ac8732SStefano Zampini     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
14459ac8732SStefano Zampini   } else { /* either full or lower-triangular (not packed) */
14559ac8732SStefano Zampini     char ord[2];
14659ac8732SStefano Zampini     if (mumps->id.ICNTL(19) == 2 || mumps->id.ICNTL(19) == 3) { /* lower triangular stored by columns or full matrix */
14759ac8732SStefano Zampini       sprintf(ord,"L");
14859ac8732SStefano Zampini     } else { /* ICNTL(19) == 1 lower triangular stored by rows */
14959ac8732SStefano Zampini       sprintf(ord,"U");
15059ac8732SStefano Zampini     }
151a0b0af32SStefano Zampini     if (mumps->schur_sym == 2) {
15259ac8732SStefano Zampini       if (!mumps->schur_pivots) {
15359ac8732SStefano Zampini         PetscScalar  lwork;
15459ac8732SStefano Zampini 
15559ac8732SStefano Zampini         ierr = PetscMalloc1(B_N,&mumps->schur_pivots);CHKERRQ(ierr);
15659ac8732SStefano Zampini         mumps->schur_B_lwork=-1;
15759ac8732SStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
15859ac8732SStefano Zampini         PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_(ord,&B_N,(PetscScalar*)mumps->id.schur,&B_slda,mumps->schur_pivots,&lwork,&mumps->schur_B_lwork,&B_ierr));
15959ac8732SStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
16059ac8732SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYTRF Lapack routine %d",(int)B_ierr);
16159ac8732SStefano Zampini         ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&mumps->schur_B_lwork);CHKERRQ(ierr);
16259ac8732SStefano Zampini         ierr = PetscMalloc1(mumps->schur_B_lwork,&mumps->schur_work);CHKERRQ(ierr);
16359ac8732SStefano Zampini       }
16459ac8732SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
16559ac8732SStefano 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));
16659ac8732SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
16759ac8732SStefano Zampini       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
16859ac8732SStefano Zampini     } else {
16959ac8732SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
17059ac8732SStefano Zampini       PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_(ord,&B_N,(PetscScalar*)mumps->id.schur,&B_slda,&B_ierr));
17159ac8732SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
17259ac8732SStefano Zampini       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
17359ac8732SStefano Zampini     }
17459ac8732SStefano Zampini   }
17559ac8732SStefano Zampini   mumps->schur_factored = PETSC_TRUE;
17659ac8732SStefano Zampini   PetscFunctionReturn(0);
17759ac8732SStefano Zampini }
17859ac8732SStefano Zampini 
17959ac8732SStefano Zampini #undef __FUNCT__
18059ac8732SStefano Zampini #define __FUNCT__ "MatMumpsInvertSchur_Private"
18159ac8732SStefano Zampini static PetscErrorCode MatMumpsInvertSchur_Private(Mat_MUMPS* mumps)
18259ac8732SStefano Zampini {
18359ac8732SStefano Zampini   PetscBLASInt   B_N,B_ierr,B_slda;
18459ac8732SStefano Zampini   PetscErrorCode ierr;
18559ac8732SStefano Zampini 
18659ac8732SStefano Zampini   PetscFunctionBegin;
18759ac8732SStefano Zampini   ierr = MatMumpsFactorSchur_Private(mumps);CHKERRQ(ierr);
18859ac8732SStefano Zampini   ierr = PetscBLASIntCast(mumps->id.size_schur,&B_N);CHKERRQ(ierr);
18959ac8732SStefano Zampini   ierr = PetscBLASIntCast(mumps->id.schur_lld,&B_slda);CHKERRQ(ierr);
19059ac8732SStefano Zampini   if (!mumps->sym) { /* MUMPS always return a full Schur matrix */
19159ac8732SStefano Zampini     if (!mumps->schur_work) {
19259ac8732SStefano Zampini       PetscScalar lwork;
19359ac8732SStefano Zampini 
19459ac8732SStefano Zampini       mumps->schur_B_lwork = -1;
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,&lwork,&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 query to GETRI Lapack routine %d",(int)B_ierr);
19959ac8732SStefano Zampini       ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&mumps->schur_B_lwork);CHKERRQ(ierr);
20059ac8732SStefano Zampini       ierr = PetscMalloc1(mumps->schur_B_lwork,&mumps->schur_work);CHKERRQ(ierr);
20159ac8732SStefano Zampini     }
20259ac8732SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
20359ac8732SStefano Zampini     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,(PetscScalar*)mumps->id.schur,&B_N,mumps->schur_pivots,mumps->schur_work,&mumps->schur_B_lwork,&B_ierr));
20459ac8732SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
20559ac8732SStefano Zampini     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
20659ac8732SStefano Zampini   } else { /* either full or lower-triangular (not packed) */
20759ac8732SStefano Zampini     char ord[2];
20859ac8732SStefano Zampini     if (mumps->id.ICNTL(19) == 2 || mumps->id.ICNTL(19) == 3) { /* lower triangular stored by columns or full matrix */
20959ac8732SStefano Zampini       sprintf(ord,"L");
21059ac8732SStefano Zampini     } else { /* ICNTL(19) == 1 lower triangular stored by rows */
21159ac8732SStefano Zampini       sprintf(ord,"U");
21259ac8732SStefano Zampini     }
213a0b0af32SStefano Zampini     if (mumps->schur_sym == 2) {
21459ac8732SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
21559ac8732SStefano Zampini       PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_(ord,&B_N,(PetscScalar*)mumps->id.schur,&B_N,mumps->schur_pivots,mumps->schur_work,&B_ierr));
21659ac8732SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
21759ac8732SStefano Zampini       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
21859ac8732SStefano Zampini     } else {
21959ac8732SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
22059ac8732SStefano Zampini       PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_(ord,&B_N,(PetscScalar*)mumps->id.schur,&B_N,&B_ierr));
22159ac8732SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
22259ac8732SStefano Zampini       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
22359ac8732SStefano Zampini     }
22459ac8732SStefano Zampini   }
22559ac8732SStefano Zampini   mumps->schur_inverted = PETSC_TRUE;
22659ac8732SStefano Zampini   PetscFunctionReturn(0);
22759ac8732SStefano Zampini }
22859ac8732SStefano Zampini 
22959ac8732SStefano Zampini #undef __FUNCT__
23059ac8732SStefano Zampini #define __FUNCT__ "MatMumpsSolveSchur_Private"
231e807eca7SStefano Zampini static PetscErrorCode MatMumpsSolveSchur_Private(Mat_MUMPS* mumps, PetscBool sol_in_redrhs)
23259ac8732SStefano Zampini {
23359ac8732SStefano Zampini   PetscBLASInt   B_N,B_Nrhs,B_ierr,B_slda,B_rlda;
23459ac8732SStefano Zampini   PetscScalar    one=1.,zero=0.;
23559ac8732SStefano Zampini   PetscErrorCode ierr;
23659ac8732SStefano Zampini 
23759ac8732SStefano Zampini   PetscFunctionBegin;
23859ac8732SStefano Zampini   ierr = MatMumpsFactorSchur_Private(mumps);CHKERRQ(ierr);
239b5fa320bSStefano Zampini   ierr = PetscBLASIntCast(mumps->id.size_schur,&B_N);CHKERRQ(ierr);
240b5fa320bSStefano Zampini   ierr = PetscBLASIntCast(mumps->id.schur_lld,&B_slda);CHKERRQ(ierr);
241b5fa320bSStefano Zampini   ierr = PetscBLASIntCast(mumps->id.nrhs,&B_Nrhs);CHKERRQ(ierr);
242b5fa320bSStefano Zampini   ierr = PetscBLASIntCast(mumps->id.lredrhs,&B_rlda);CHKERRQ(ierr);
24359ac8732SStefano Zampini   if (mumps->schur_inverted) {
24459ac8732SStefano Zampini     PetscInt sizesol = B_Nrhs*B_N;
24559ac8732SStefano Zampini     if (!mumps->schur_sol || sizesol > mumps->schur_sizesol) {
24659ac8732SStefano Zampini       ierr = PetscFree(mumps->schur_sol);CHKERRQ(ierr);
24759ac8732SStefano Zampini       ierr = PetscMalloc1(sizesol,&mumps->schur_sol);CHKERRQ(ierr);
24859ac8732SStefano Zampini       mumps->schur_sizesol = sizesol;
249b5fa320bSStefano Zampini     }
25059ac8732SStefano Zampini     if (!mumps->sym) {
25159ac8732SStefano Zampini       char type[2];
252b5fa320bSStefano Zampini       if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
25359ac8732SStefano Zampini         if (!mumps->id.ICNTL(9)) { /* transpose solve */
25459ac8732SStefano Zampini           sprintf(type,"N");
255b5fa320bSStefano Zampini         } else {
25659ac8732SStefano Zampini           sprintf(type,"T");
257b5fa320bSStefano Zampini         }
25859ac8732SStefano Zampini       } else { /* stored by columns */
25959ac8732SStefano Zampini         if (!mumps->id.ICNTL(9)) { /* transpose solve */
26059ac8732SStefano Zampini           sprintf(type,"T");
26159ac8732SStefano Zampini         } else {
26259ac8732SStefano Zampini           sprintf(type,"N");
26359ac8732SStefano Zampini         }
26459ac8732SStefano Zampini       }
26559ac8732SStefano 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));
26659ac8732SStefano Zampini     } else {
26759ac8732SStefano Zampini       char ord[2];
26859ac8732SStefano Zampini       if (mumps->id.ICNTL(19) == 2 || mumps->id.ICNTL(19) == 3) { /* lower triangular stored by columns or full matrix */
26959ac8732SStefano Zampini         sprintf(ord,"L");
27059ac8732SStefano Zampini       } else { /* ICNTL(19) == 1 lower triangular stored by rows */
27159ac8732SStefano Zampini         sprintf(ord,"U");
27259ac8732SStefano Zampini       }
27359ac8732SStefano 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));
27459ac8732SStefano Zampini     }
275e807eca7SStefano Zampini     if (sol_in_redrhs) {
27659ac8732SStefano Zampini       ierr = PetscMemcpy(mumps->id.redrhs,mumps->schur_sol,sizesol*sizeof(PetscScalar));CHKERRQ(ierr);
277e807eca7SStefano Zampini     }
2783280fbe9SStefano Zampini   } else { /* Schur complement has not yet been inverted */
279a12f35bfSStefano Zampini     MumpsScalar *orhs=NULL;
280a12f35bfSStefano Zampini 
281a12f35bfSStefano Zampini     if (!sol_in_redrhs) {
282a12f35bfSStefano Zampini       PetscInt sizesol = B_Nrhs*B_N;
283a12f35bfSStefano Zampini       if (!mumps->schur_sol || sizesol > mumps->schur_sizesol) {
284a12f35bfSStefano Zampini         ierr = PetscFree(mumps->schur_sol);CHKERRQ(ierr);
285a12f35bfSStefano Zampini         ierr = PetscMalloc1(sizesol,&mumps->schur_sol);CHKERRQ(ierr);
286a12f35bfSStefano Zampini         mumps->schur_sizesol = sizesol;
287a12f35bfSStefano Zampini       }
288a12f35bfSStefano Zampini       orhs = mumps->id.redrhs;
289a12f35bfSStefano Zampini       ierr = PetscMemcpy(mumps->schur_sol,mumps->id.redrhs,sizesol*sizeof(PetscScalar));CHKERRQ(ierr);
290a12f35bfSStefano Zampini       mumps->id.redrhs = (MumpsScalar*)mumps->schur_sol;
291a12f35bfSStefano Zampini     }
29259ac8732SStefano Zampini     if (!mumps->sym) { /* MUMPS always return a full Schur matrix */
29359ac8732SStefano Zampini       char type[2];
29459ac8732SStefano Zampini       if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
29559ac8732SStefano Zampini         if (!mumps->id.ICNTL(9)) { /* transpose solve */
29659ac8732SStefano Zampini           sprintf(type,"N");
29759ac8732SStefano Zampini         } else {
29859ac8732SStefano Zampini           sprintf(type,"T");
29959ac8732SStefano Zampini         }
30059ac8732SStefano Zampini       } else { /* stored by columns */
30159ac8732SStefano Zampini         if (!mumps->id.ICNTL(9)) { /* transpose solve */
30259ac8732SStefano Zampini           sprintf(type,"T");
30359ac8732SStefano Zampini         } else {
30459ac8732SStefano Zampini           sprintf(type,"N");
30559ac8732SStefano Zampini         }
30659ac8732SStefano Zampini       }
30759ac8732SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
30859ac8732SStefano 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));
309b5fa320bSStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
310b5fa320bSStefano Zampini       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRS Lapack routine %d",(int)B_ierr);
311b5fa320bSStefano Zampini     } else { /* either full or lower-triangular (not packed) */
312b5fa320bSStefano Zampini       char ord[2];
313b5fa320bSStefano Zampini       if (mumps->id.ICNTL(19) == 2 || mumps->id.ICNTL(19) == 3) { /* lower triangular stored by columns or full matrix */
314b5fa320bSStefano Zampini         sprintf(ord,"L");
315b5fa320bSStefano Zampini       } else { /* ICNTL(19) == 1 lower triangular stored by rows */
316b5fa320bSStefano Zampini         sprintf(ord,"U");
317b5fa320bSStefano Zampini       }
318a0b0af32SStefano Zampini       if (mumps->schur_sym == 2) {
319b5fa320bSStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
32059ac8732SStefano 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));
321b5fa320bSStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
322b5fa320bSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRS Lapack routine %d",(int)B_ierr);
323b5fa320bSStefano Zampini       } else {
324b5fa320bSStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
32559ac8732SStefano Zampini         PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_(ord,&B_N,&B_Nrhs,(PetscScalar*)mumps->id.schur,&B_slda,(PetscScalar*)mumps->id.redrhs,&B_rlda,&B_ierr));
326b5fa320bSStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
327b5fa320bSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRS Lapack routine %d",(int)B_ierr);
328b5fa320bSStefano Zampini       }
329b5fa320bSStefano Zampini     }
330e807eca7SStefano Zampini     if (!sol_in_redrhs) {
331a12f35bfSStefano Zampini       mumps->id.redrhs = orhs;
332e807eca7SStefano Zampini     }
33359ac8732SStefano Zampini   }
334b5fa320bSStefano Zampini   PetscFunctionReturn(0);
335b5fa320bSStefano Zampini }
336b5fa320bSStefano Zampini 
33759ac8732SStefano Zampini #undef __FUNCT__
33859ac8732SStefano Zampini #define __FUNCT__ "MatMumpsHandleSchur_Private"
339b8f61ee1SStefano Zampini static PetscErrorCode MatMumpsHandleSchur_Private(Mat_MUMPS* mumps, PetscBool expansion)
340b5fa320bSStefano Zampini {
341b5fa320bSStefano Zampini   PetscErrorCode ierr;
342b5fa320bSStefano Zampini 
343b5fa320bSStefano Zampini   PetscFunctionBegin;
344b5fa320bSStefano Zampini   if (!mumps->id.ICNTL(19)) { /* do nothing when Schur complement has not been computed */
345b5fa320bSStefano Zampini     PetscFunctionReturn(0);
346b5fa320bSStefano Zampini   }
347b8f61ee1SStefano Zampini   if (!expansion) { /* prepare for the condensation step */
348b5fa320bSStefano Zampini     PetscInt sizeredrhs = mumps->id.nrhs*mumps->id.size_schur;
349b5fa320bSStefano Zampini     /* allocate MUMPS internal array to store reduced right-hand sides */
350b5fa320bSStefano Zampini     if (!mumps->id.redrhs || sizeredrhs > mumps->sizeredrhs) {
351b5fa320bSStefano Zampini       ierr = PetscFree(mumps->id.redrhs);CHKERRQ(ierr);
352b5fa320bSStefano Zampini       mumps->id.lredrhs = mumps->id.size_schur;
353b5fa320bSStefano Zampini       ierr = PetscMalloc1(mumps->id.nrhs*mumps->id.lredrhs,&mumps->id.redrhs);CHKERRQ(ierr);
354b5fa320bSStefano Zampini       mumps->sizeredrhs = mumps->id.nrhs*mumps->id.lredrhs;
355b5fa320bSStefano Zampini     }
356b5fa320bSStefano Zampini     mumps->id.ICNTL(26) = 1; /* condensation phase */
357b5fa320bSStefano Zampini   } else { /* prepare for the expansion step */
358b8f61ee1SStefano Zampini     /* solve Schur complement (this has to be done by the MUMPS user, so basically us) */
359e807eca7SStefano Zampini     ierr = MatMumpsSolveSchur_Private(mumps,PETSC_TRUE);CHKERRQ(ierr);
360b5fa320bSStefano Zampini     mumps->id.ICNTL(26) = 2; /* expansion phase */
361b5fa320bSStefano Zampini     PetscMUMPS_c(&mumps->id);
362b5fa320bSStefano 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));
363b5fa320bSStefano Zampini     /* restore defaults */
364b5fa320bSStefano Zampini     mumps->id.ICNTL(26) = -1;
365d3d598ffSStefano Zampini     /* free MUMPS internal array for redrhs if we have solved for multiple rhs in order to save memory space */
366d3d598ffSStefano Zampini     if (mumps->id.nrhs > 1) {
367d3d598ffSStefano Zampini       ierr = PetscFree(mumps->id.redrhs);CHKERRQ(ierr);
368d3d598ffSStefano Zampini       mumps->id.lredrhs = 0;
369d3d598ffSStefano Zampini       mumps->sizeredrhs = 0;
370d3d598ffSStefano Zampini     }
371b5fa320bSStefano Zampini   }
372b5fa320bSStefano Zampini   PetscFunctionReturn(0);
373b5fa320bSStefano Zampini }
374b5fa320bSStefano Zampini 
375397b6df1SKris Buschelman /*
376d341cd04SHong Zhang   MatConvertToTriples_A_B - convert Petsc matrix to triples: row[nz], col[nz], val[nz]
377d341cd04SHong Zhang 
378397b6df1SKris Buschelman   input:
37967877ebaSShri Abhyankar     A       - matrix in aij,baij or sbaij (bs=1) format
380397b6df1SKris Buschelman     shift   - 0: C style output triple; 1: Fortran style output triple.
381bccb9932SShri Abhyankar     reuse   - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple
382bccb9932SShri Abhyankar               MAT_REUSE_MATRIX:   only the values in v array are updated
383397b6df1SKris Buschelman   output:
384397b6df1SKris Buschelman     nnz     - dim of r, c, and v (number of local nonzero entries of A)
385397b6df1SKris Buschelman     r, c, v - row and col index, matrix values (matrix triples)
386eb9baa12SBarry Smith 
387eb9baa12SBarry Smith   The returned values r, c, and sometimes v are obtained in a single PetscMalloc(). Then in MatDestroy_MUMPS() it is
388eb9baa12SBarry Smith   freed with PetscFree((mumps->irn);  This is not ideal code, the fact that v is ONLY sometimes part of mumps->irn means
389eb9baa12SBarry Smith   that the PetscMalloc() cannot easily be replaced with a PetscMalloc3().
390eb9baa12SBarry Smith 
391397b6df1SKris Buschelman  */
39216ebf90aSShri Abhyankar 
39316ebf90aSShri Abhyankar #undef __FUNCT__
39416ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqaij"
395bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
396b24902e0SBarry Smith {
397185f6596SHong Zhang   const PetscInt *ai,*aj,*ajj,M=A->rmap->n;
39867877ebaSShri Abhyankar   PetscInt       nz,rnz,i,j;
399dfbe8321SBarry Smith   PetscErrorCode ierr;
400c1490034SHong Zhang   PetscInt       *row,*col;
40116ebf90aSShri Abhyankar   Mat_SeqAIJ     *aa=(Mat_SeqAIJ*)A->data;
402397b6df1SKris Buschelman 
403397b6df1SKris Buschelman   PetscFunctionBegin;
40416ebf90aSShri Abhyankar   *v=aa->a;
405bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
4062205254eSKarl Rupp     nz   = aa->nz;
4072205254eSKarl Rupp     ai   = aa->i;
4082205254eSKarl Rupp     aj   = aa->j;
40916ebf90aSShri Abhyankar     *nnz = nz;
410785e854fSJed Brown     ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr);
411185f6596SHong Zhang     col  = row + nz;
412185f6596SHong Zhang 
41316ebf90aSShri Abhyankar     nz = 0;
41416ebf90aSShri Abhyankar     for (i=0; i<M; i++) {
41516ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
41667877ebaSShri Abhyankar       ajj = aj + ai[i];
41767877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
41867877ebaSShri Abhyankar         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
41916ebf90aSShri Abhyankar       }
42016ebf90aSShri Abhyankar     }
42116ebf90aSShri Abhyankar     *r = row; *c = col;
42216ebf90aSShri Abhyankar   }
42316ebf90aSShri Abhyankar   PetscFunctionReturn(0);
42416ebf90aSShri Abhyankar }
425397b6df1SKris Buschelman 
42616ebf90aSShri Abhyankar #undef __FUNCT__
42767877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqbaij_seqaij"
428bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
42967877ebaSShri Abhyankar {
43067877ebaSShri Abhyankar   Mat_SeqBAIJ    *aa=(Mat_SeqBAIJ*)A->data;
43133d57670SJed Brown   const PetscInt *ai,*aj,*ajj,bs2 = aa->bs2;
43233d57670SJed Brown   PetscInt       bs,M,nz,idx=0,rnz,i,j,k,m;
43367877ebaSShri Abhyankar   PetscErrorCode ierr;
43467877ebaSShri Abhyankar   PetscInt       *row,*col;
43567877ebaSShri Abhyankar 
43667877ebaSShri Abhyankar   PetscFunctionBegin;
43733d57670SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
43833d57670SJed Brown   M = A->rmap->N/bs;
439cf3759fdSShri Abhyankar   *v = aa->a;
440bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
441cf3759fdSShri Abhyankar     ai   = aa->i; aj = aa->j;
44267877ebaSShri Abhyankar     nz   = bs2*aa->nz;
44367877ebaSShri Abhyankar     *nnz = nz;
444785e854fSJed Brown     ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr);
445185f6596SHong Zhang     col  = row + nz;
446185f6596SHong Zhang 
44767877ebaSShri Abhyankar     for (i=0; i<M; i++) {
44867877ebaSShri Abhyankar       ajj = aj + ai[i];
44967877ebaSShri Abhyankar       rnz = ai[i+1] - ai[i];
45067877ebaSShri Abhyankar       for (k=0; k<rnz; k++) {
45167877ebaSShri Abhyankar         for (j=0; j<bs; j++) {
45267877ebaSShri Abhyankar           for (m=0; m<bs; m++) {
45367877ebaSShri Abhyankar             row[idx]   = i*bs + m + shift;
454cf3759fdSShri Abhyankar             col[idx++] = bs*(ajj[k]) + j + shift;
45567877ebaSShri Abhyankar           }
45667877ebaSShri Abhyankar         }
45767877ebaSShri Abhyankar       }
45867877ebaSShri Abhyankar     }
459cf3759fdSShri Abhyankar     *r = row; *c = col;
46067877ebaSShri Abhyankar   }
46167877ebaSShri Abhyankar   PetscFunctionReturn(0);
46267877ebaSShri Abhyankar }
46367877ebaSShri Abhyankar 
46467877ebaSShri Abhyankar #undef __FUNCT__
46516ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqsbaij_seqsbaij"
466bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
46716ebf90aSShri Abhyankar {
46867877ebaSShri Abhyankar   const PetscInt *ai, *aj,*ajj,M=A->rmap->n;
46967877ebaSShri Abhyankar   PetscInt       nz,rnz,i,j;
47016ebf90aSShri Abhyankar   PetscErrorCode ierr;
47116ebf90aSShri Abhyankar   PetscInt       *row,*col;
47216ebf90aSShri Abhyankar   Mat_SeqSBAIJ   *aa=(Mat_SeqSBAIJ*)A->data;
47316ebf90aSShri Abhyankar 
47416ebf90aSShri Abhyankar   PetscFunctionBegin;
475882afa5aSHong Zhang   *v = aa->a;
476bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
4772205254eSKarl Rupp     nz   = aa->nz;
4782205254eSKarl Rupp     ai   = aa->i;
4792205254eSKarl Rupp     aj   = aa->j;
4802205254eSKarl Rupp     *v   = aa->a;
48116ebf90aSShri Abhyankar     *nnz = nz;
482785e854fSJed Brown     ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr);
483185f6596SHong Zhang     col  = row + nz;
484185f6596SHong Zhang 
48516ebf90aSShri Abhyankar     nz = 0;
48616ebf90aSShri Abhyankar     for (i=0; i<M; i++) {
48716ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
48867877ebaSShri Abhyankar       ajj = aj + ai[i];
48967877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
49067877ebaSShri Abhyankar         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
49116ebf90aSShri Abhyankar       }
49216ebf90aSShri Abhyankar     }
49316ebf90aSShri Abhyankar     *r = row; *c = col;
49416ebf90aSShri Abhyankar   }
49516ebf90aSShri Abhyankar   PetscFunctionReturn(0);
49616ebf90aSShri Abhyankar }
49716ebf90aSShri Abhyankar 
49816ebf90aSShri Abhyankar #undef __FUNCT__
49916ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqsbaij"
500bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
50116ebf90aSShri Abhyankar {
50267877ebaSShri Abhyankar   const PetscInt    *ai,*aj,*ajj,*adiag,M=A->rmap->n;
50367877ebaSShri Abhyankar   PetscInt          nz,rnz,i,j;
50467877ebaSShri Abhyankar   const PetscScalar *av,*v1;
50516ebf90aSShri Abhyankar   PetscScalar       *val;
50616ebf90aSShri Abhyankar   PetscErrorCode    ierr;
50716ebf90aSShri Abhyankar   PetscInt          *row,*col;
508829b1710SHong Zhang   Mat_SeqAIJ        *aa=(Mat_SeqAIJ*)A->data;
50916ebf90aSShri Abhyankar 
51016ebf90aSShri Abhyankar   PetscFunctionBegin;
51116ebf90aSShri Abhyankar   ai   =aa->i; aj=aa->j;av=aa->a;
51216ebf90aSShri Abhyankar   adiag=aa->diag;
513bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
514829b1710SHong Zhang     /* count nz in the uppper triangular part of A */
515829b1710SHong Zhang     nz = 0;
516829b1710SHong Zhang     for (i=0; i<M; i++) nz += ai[i+1] - adiag[i];
51716ebf90aSShri Abhyankar     *nnz = nz;
518829b1710SHong Zhang 
519185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
520185f6596SHong Zhang     col  = row + nz;
521185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
522185f6596SHong Zhang 
52316ebf90aSShri Abhyankar     nz = 0;
52416ebf90aSShri Abhyankar     for (i=0; i<M; i++) {
52516ebf90aSShri Abhyankar       rnz = ai[i+1] - adiag[i];
52667877ebaSShri Abhyankar       ajj = aj + adiag[i];
527cf3759fdSShri Abhyankar       v1  = av + adiag[i];
52867877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
52967877ebaSShri Abhyankar         row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
53016ebf90aSShri Abhyankar       }
53116ebf90aSShri Abhyankar     }
53216ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
533397b6df1SKris Buschelman   } else {
53416ebf90aSShri Abhyankar     nz = 0; val = *v;
53516ebf90aSShri Abhyankar     for (i=0; i <M; i++) {
53616ebf90aSShri Abhyankar       rnz = ai[i+1] - adiag[i];
53767877ebaSShri Abhyankar       ajj = aj + adiag[i];
53867877ebaSShri Abhyankar       v1  = av + adiag[i];
53967877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
54067877ebaSShri Abhyankar         val[nz++] = v1[j];
54116ebf90aSShri Abhyankar       }
54216ebf90aSShri Abhyankar     }
54316ebf90aSShri Abhyankar   }
54416ebf90aSShri Abhyankar   PetscFunctionReturn(0);
54516ebf90aSShri Abhyankar }
54616ebf90aSShri Abhyankar 
54716ebf90aSShri Abhyankar #undef __FUNCT__
54816ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpisbaij_mpisbaij"
549bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
55016ebf90aSShri Abhyankar {
55116ebf90aSShri Abhyankar   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
55216ebf90aSShri Abhyankar   PetscErrorCode    ierr;
55316ebf90aSShri Abhyankar   PetscInt          rstart,nz,i,j,jj,irow,countA,countB;
55416ebf90aSShri Abhyankar   PetscInt          *row,*col;
55516ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
55616ebf90aSShri Abhyankar   PetscScalar       *val;
557397b6df1SKris Buschelman   Mat_MPISBAIJ      *mat = (Mat_MPISBAIJ*)A->data;
558397b6df1SKris Buschelman   Mat_SeqSBAIJ      *aa  = (Mat_SeqSBAIJ*)(mat->A)->data;
559397b6df1SKris Buschelman   Mat_SeqBAIJ       *bb  = (Mat_SeqBAIJ*)(mat->B)->data;
56016ebf90aSShri Abhyankar 
56116ebf90aSShri Abhyankar   PetscFunctionBegin;
562d0f46423SBarry Smith   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
563397b6df1SKris Buschelman   av=aa->a; bv=bb->a;
564397b6df1SKris Buschelman 
5652205254eSKarl Rupp   garray = mat->garray;
5662205254eSKarl Rupp 
567bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
56816ebf90aSShri Abhyankar     nz   = aa->nz + bb->nz;
56916ebf90aSShri Abhyankar     *nnz = nz;
570185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
571185f6596SHong Zhang     col  = row + nz;
572185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
573185f6596SHong Zhang 
574397b6df1SKris Buschelman     *r = row; *c = col; *v = val;
575397b6df1SKris Buschelman   } else {
576397b6df1SKris Buschelman     row = *r; col = *c; val = *v;
577397b6df1SKris Buschelman   }
578397b6df1SKris Buschelman 
579028e57e8SHong Zhang   jj = 0; irow = rstart;
580397b6df1SKris Buschelman   for (i=0; i<m; i++) {
581397b6df1SKris Buschelman     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
582397b6df1SKris Buschelman     countA = ai[i+1] - ai[i];
583397b6df1SKris Buschelman     countB = bi[i+1] - bi[i];
584397b6df1SKris Buschelman     bjj    = bj + bi[i];
58516ebf90aSShri Abhyankar     v1     = av + ai[i];
58616ebf90aSShri Abhyankar     v2     = bv + bi[i];
587397b6df1SKris Buschelman 
588397b6df1SKris Buschelman     /* A-part */
589397b6df1SKris Buschelman     for (j=0; j<countA; j++) {
590bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
591397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
592397b6df1SKris Buschelman       }
59316ebf90aSShri Abhyankar       val[jj++] = v1[j];
594397b6df1SKris Buschelman     }
59516ebf90aSShri Abhyankar 
59616ebf90aSShri Abhyankar     /* B-part */
59716ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
598bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
599397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
600397b6df1SKris Buschelman       }
60116ebf90aSShri Abhyankar       val[jj++] = v2[j];
60216ebf90aSShri Abhyankar     }
60316ebf90aSShri Abhyankar     irow++;
60416ebf90aSShri Abhyankar   }
60516ebf90aSShri Abhyankar   PetscFunctionReturn(0);
60616ebf90aSShri Abhyankar }
60716ebf90aSShri Abhyankar 
60816ebf90aSShri Abhyankar #undef __FUNCT__
60916ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpiaij"
610bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
61116ebf90aSShri Abhyankar {
61216ebf90aSShri Abhyankar   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
61316ebf90aSShri Abhyankar   PetscErrorCode    ierr;
61416ebf90aSShri Abhyankar   PetscInt          rstart,nz,i,j,jj,irow,countA,countB;
61516ebf90aSShri Abhyankar   PetscInt          *row,*col;
61616ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
61716ebf90aSShri Abhyankar   PetscScalar       *val;
61816ebf90aSShri Abhyankar   Mat_MPIAIJ        *mat = (Mat_MPIAIJ*)A->data;
61916ebf90aSShri Abhyankar   Mat_SeqAIJ        *aa  = (Mat_SeqAIJ*)(mat->A)->data;
62016ebf90aSShri Abhyankar   Mat_SeqAIJ        *bb  = (Mat_SeqAIJ*)(mat->B)->data;
62116ebf90aSShri Abhyankar 
62216ebf90aSShri Abhyankar   PetscFunctionBegin;
62316ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
62416ebf90aSShri Abhyankar   av=aa->a; bv=bb->a;
62516ebf90aSShri Abhyankar 
6262205254eSKarl Rupp   garray = mat->garray;
6272205254eSKarl Rupp 
628bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
62916ebf90aSShri Abhyankar     nz   = aa->nz + bb->nz;
63016ebf90aSShri Abhyankar     *nnz = nz;
631185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
632185f6596SHong Zhang     col  = row + nz;
633185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
634185f6596SHong Zhang 
63516ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
63616ebf90aSShri Abhyankar   } else {
63716ebf90aSShri Abhyankar     row = *r; col = *c; val = *v;
63816ebf90aSShri Abhyankar   }
63916ebf90aSShri Abhyankar 
64016ebf90aSShri Abhyankar   jj = 0; irow = rstart;
64116ebf90aSShri Abhyankar   for (i=0; i<m; i++) {
64216ebf90aSShri Abhyankar     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
64316ebf90aSShri Abhyankar     countA = ai[i+1] - ai[i];
64416ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
64516ebf90aSShri Abhyankar     bjj    = bj + bi[i];
64616ebf90aSShri Abhyankar     v1     = av + ai[i];
64716ebf90aSShri Abhyankar     v2     = bv + bi[i];
64816ebf90aSShri Abhyankar 
64916ebf90aSShri Abhyankar     /* A-part */
65016ebf90aSShri Abhyankar     for (j=0; j<countA; j++) {
651bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
65216ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
65316ebf90aSShri Abhyankar       }
65416ebf90aSShri Abhyankar       val[jj++] = v1[j];
65516ebf90aSShri Abhyankar     }
65616ebf90aSShri Abhyankar 
65716ebf90aSShri Abhyankar     /* B-part */
65816ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
659bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
66016ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
66116ebf90aSShri Abhyankar       }
66216ebf90aSShri Abhyankar       val[jj++] = v2[j];
66316ebf90aSShri Abhyankar     }
66416ebf90aSShri Abhyankar     irow++;
66516ebf90aSShri Abhyankar   }
66616ebf90aSShri Abhyankar   PetscFunctionReturn(0);
66716ebf90aSShri Abhyankar }
66816ebf90aSShri Abhyankar 
66916ebf90aSShri Abhyankar #undef __FUNCT__
67067877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpibaij_mpiaij"
671bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
67267877ebaSShri Abhyankar {
67367877ebaSShri Abhyankar   Mat_MPIBAIJ       *mat    = (Mat_MPIBAIJ*)A->data;
67467877ebaSShri Abhyankar   Mat_SeqBAIJ       *aa     = (Mat_SeqBAIJ*)(mat->A)->data;
67567877ebaSShri Abhyankar   Mat_SeqBAIJ       *bb     = (Mat_SeqBAIJ*)(mat->B)->data;
67667877ebaSShri Abhyankar   const PetscInt    *ai     = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj;
677d985c460SShri Abhyankar   const PetscInt    *garray = mat->garray,mbs=mat->mbs,rstart=A->rmap->rstart;
67833d57670SJed Brown   const PetscInt    bs2=mat->bs2;
67967877ebaSShri Abhyankar   PetscErrorCode    ierr;
68033d57670SJed Brown   PetscInt          bs,nz,i,j,k,n,jj,irow,countA,countB,idx;
68167877ebaSShri Abhyankar   PetscInt          *row,*col;
68267877ebaSShri Abhyankar   const PetscScalar *av=aa->a, *bv=bb->a,*v1,*v2;
68367877ebaSShri Abhyankar   PetscScalar       *val;
68467877ebaSShri Abhyankar 
68567877ebaSShri Abhyankar   PetscFunctionBegin;
68633d57670SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
687bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
68867877ebaSShri Abhyankar     nz   = bs2*(aa->nz + bb->nz);
68967877ebaSShri Abhyankar     *nnz = nz;
690185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
691185f6596SHong Zhang     col  = row + nz;
692185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
693185f6596SHong Zhang 
69467877ebaSShri Abhyankar     *r = row; *c = col; *v = val;
69567877ebaSShri Abhyankar   } else {
69667877ebaSShri Abhyankar     row = *r; col = *c; val = *v;
69767877ebaSShri Abhyankar   }
69867877ebaSShri Abhyankar 
699d985c460SShri Abhyankar   jj = 0; irow = rstart;
70067877ebaSShri Abhyankar   for (i=0; i<mbs; i++) {
70167877ebaSShri Abhyankar     countA = ai[i+1] - ai[i];
70267877ebaSShri Abhyankar     countB = bi[i+1] - bi[i];
70367877ebaSShri Abhyankar     ajj    = aj + ai[i];
70467877ebaSShri Abhyankar     bjj    = bj + bi[i];
70567877ebaSShri Abhyankar     v1     = av + bs2*ai[i];
70667877ebaSShri Abhyankar     v2     = bv + bs2*bi[i];
70767877ebaSShri Abhyankar 
70867877ebaSShri Abhyankar     idx = 0;
70967877ebaSShri Abhyankar     /* A-part */
71067877ebaSShri Abhyankar     for (k=0; k<countA; k++) {
71167877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
71267877ebaSShri Abhyankar         for (n=0; n<bs; n++) {
713bccb9932SShri Abhyankar           if (reuse == MAT_INITIAL_MATRIX) {
714d985c460SShri Abhyankar             row[jj] = irow + n + shift;
715d985c460SShri Abhyankar             col[jj] = rstart + bs*ajj[k] + j + shift;
71667877ebaSShri Abhyankar           }
71767877ebaSShri Abhyankar           val[jj++] = v1[idx++];
71867877ebaSShri Abhyankar         }
71967877ebaSShri Abhyankar       }
72067877ebaSShri Abhyankar     }
72167877ebaSShri Abhyankar 
72267877ebaSShri Abhyankar     idx = 0;
72367877ebaSShri Abhyankar     /* B-part */
72467877ebaSShri Abhyankar     for (k=0; k<countB; k++) {
72567877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
72667877ebaSShri Abhyankar         for (n=0; n<bs; n++) {
727bccb9932SShri Abhyankar           if (reuse == MAT_INITIAL_MATRIX) {
728d985c460SShri Abhyankar             row[jj] = irow + n + shift;
729d985c460SShri Abhyankar             col[jj] = bs*garray[bjj[k]] + j + shift;
73067877ebaSShri Abhyankar           }
731d985c460SShri Abhyankar           val[jj++] = v2[idx++];
73267877ebaSShri Abhyankar         }
73367877ebaSShri Abhyankar       }
73467877ebaSShri Abhyankar     }
735d985c460SShri Abhyankar     irow += bs;
73667877ebaSShri Abhyankar   }
73767877ebaSShri Abhyankar   PetscFunctionReturn(0);
73867877ebaSShri Abhyankar }
73967877ebaSShri Abhyankar 
74067877ebaSShri Abhyankar #undef __FUNCT__
74116ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpisbaij"
742bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
74316ebf90aSShri Abhyankar {
74416ebf90aSShri Abhyankar   const PetscInt    *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
74516ebf90aSShri Abhyankar   PetscErrorCode    ierr;
746e0bace9bSHong Zhang   PetscInt          rstart,nz,nza,nzb,i,j,jj,irow,countA,countB;
74716ebf90aSShri Abhyankar   PetscInt          *row,*col;
74816ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
74916ebf90aSShri Abhyankar   PetscScalar       *val;
75016ebf90aSShri Abhyankar   Mat_MPIAIJ        *mat =  (Mat_MPIAIJ*)A->data;
75116ebf90aSShri Abhyankar   Mat_SeqAIJ        *aa  =(Mat_SeqAIJ*)(mat->A)->data;
75216ebf90aSShri Abhyankar   Mat_SeqAIJ        *bb  =(Mat_SeqAIJ*)(mat->B)->data;
75316ebf90aSShri Abhyankar 
75416ebf90aSShri Abhyankar   PetscFunctionBegin;
75516ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j; adiag=aa->diag;
75616ebf90aSShri Abhyankar   bi=bb->i; bj=bb->j; garray = mat->garray;
75716ebf90aSShri Abhyankar   av=aa->a; bv=bb->a;
7582205254eSKarl Rupp 
75916ebf90aSShri Abhyankar   rstart = A->rmap->rstart;
76016ebf90aSShri Abhyankar 
761bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
762e0bace9bSHong Zhang     nza = 0;    /* num of upper triangular entries in mat->A, including diagonals */
763e0bace9bSHong Zhang     nzb = 0;    /* num of upper triangular entries in mat->B */
76416ebf90aSShri Abhyankar     for (i=0; i<m; i++) {
765e0bace9bSHong Zhang       nza   += (ai[i+1] - adiag[i]);
76616ebf90aSShri Abhyankar       countB = bi[i+1] - bi[i];
76716ebf90aSShri Abhyankar       bjj    = bj + bi[i];
768e0bace9bSHong Zhang       for (j=0; j<countB; j++) {
769e0bace9bSHong Zhang         if (garray[bjj[j]] > rstart) nzb++;
770e0bace9bSHong Zhang       }
771e0bace9bSHong Zhang     }
77216ebf90aSShri Abhyankar 
773e0bace9bSHong Zhang     nz   = nza + nzb; /* total nz of upper triangular part of mat */
77416ebf90aSShri Abhyankar     *nnz = nz;
775185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
776185f6596SHong Zhang     col  = row + nz;
777185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
778185f6596SHong Zhang 
77916ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
78016ebf90aSShri Abhyankar   } else {
78116ebf90aSShri Abhyankar     row = *r; col = *c; val = *v;
78216ebf90aSShri Abhyankar   }
78316ebf90aSShri Abhyankar 
78416ebf90aSShri Abhyankar   jj = 0; irow = rstart;
78516ebf90aSShri Abhyankar   for (i=0; i<m; i++) {
78616ebf90aSShri Abhyankar     ajj    = aj + adiag[i];                 /* ptr to the beginning of the diagonal of this row */
78716ebf90aSShri Abhyankar     v1     = av + adiag[i];
78816ebf90aSShri Abhyankar     countA = ai[i+1] - adiag[i];
78916ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
79016ebf90aSShri Abhyankar     bjj    = bj + bi[i];
79116ebf90aSShri Abhyankar     v2     = bv + bi[i];
79216ebf90aSShri Abhyankar 
79316ebf90aSShri Abhyankar     /* A-part */
79416ebf90aSShri Abhyankar     for (j=0; j<countA; j++) {
795bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
79616ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
79716ebf90aSShri Abhyankar       }
79816ebf90aSShri Abhyankar       val[jj++] = v1[j];
79916ebf90aSShri Abhyankar     }
80016ebf90aSShri Abhyankar 
80116ebf90aSShri Abhyankar     /* B-part */
80216ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
80316ebf90aSShri Abhyankar       if (garray[bjj[j]] > rstart) {
804bccb9932SShri Abhyankar         if (reuse == MAT_INITIAL_MATRIX) {
80516ebf90aSShri Abhyankar           row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
80616ebf90aSShri Abhyankar         }
80716ebf90aSShri Abhyankar         val[jj++] = v2[j];
80816ebf90aSShri Abhyankar       }
809397b6df1SKris Buschelman     }
810397b6df1SKris Buschelman     irow++;
811397b6df1SKris Buschelman   }
812397b6df1SKris Buschelman   PetscFunctionReturn(0);
813397b6df1SKris Buschelman }
814397b6df1SKris Buschelman 
815397b6df1SKris Buschelman #undef __FUNCT__
81620be8e61SHong Zhang #define __FUNCT__ "MatGetDiagonal_MUMPS"
81720be8e61SHong Zhang PetscErrorCode MatGetDiagonal_MUMPS(Mat A,Vec v)
81820be8e61SHong Zhang {
81920be8e61SHong Zhang   PetscFunctionBegin;
82020be8e61SHong Zhang   SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type: MUMPS factor");
82120be8e61SHong Zhang   PetscFunctionReturn(0);
82220be8e61SHong Zhang }
82320be8e61SHong Zhang 
82420be8e61SHong Zhang #undef __FUNCT__
8253924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS"
826dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A)
827dfbe8321SBarry Smith {
828a5e57a09SHong Zhang   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
829dfbe8321SBarry Smith   PetscErrorCode ierr;
830b24902e0SBarry Smith 
831397b6df1SKris Buschelman   PetscFunctionBegin;
832a5e57a09SHong Zhang   ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr);
833a5e57a09SHong Zhang   ierr = VecScatterDestroy(&mumps->scat_rhs);CHKERRQ(ierr);
834a5e57a09SHong Zhang   ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
835801fbe65SHong Zhang   ierr = VecDestroy(&mumps->b_seq);CHKERRQ(ierr);
836a5e57a09SHong Zhang   ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr);
837a5e57a09SHong Zhang   ierr = PetscFree(mumps->id.perm_in);CHKERRQ(ierr);
838a5e57a09SHong Zhang   ierr = PetscFree(mumps->irn);CHKERRQ(ierr);
839b34f08ffSHong Zhang   ierr = PetscFree(mumps->info);CHKERRQ(ierr);
84059ac8732SStefano Zampini   ierr = MatMumpsResetSchur_Private(mumps);CHKERRQ(ierr);
841a5e57a09SHong Zhang   mumps->id.job = JOB_END;
842a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
8436f3cc6f9SBarry Smith   ierr = MPI_Comm_free(&mumps->comm_mumps);CHKERRQ(ierr);
844a5e57a09SHong Zhang   if (mumps->Destroy) {
845a5e57a09SHong Zhang     ierr = (mumps->Destroy)(A);CHKERRQ(ierr);
846bf0cc555SLisandro Dalcin   }
847bf0cc555SLisandro Dalcin   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
848bf0cc555SLisandro Dalcin 
84997969023SHong Zhang   /* clear composed functions */
850bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr);
8515a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorSetSchurIS_C",NULL);CHKERRQ(ierr);
8525a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorInvertSchurComplement_C",NULL);CHKERRQ(ierr);
8535a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorCreateSchurComplement_C",NULL);CHKERRQ(ierr);
8545a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSchurComplement_C",NULL);CHKERRQ(ierr);
8555a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorSolveSchurComplement_C",NULL);CHKERRQ(ierr);
8565a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorSolveSchurComplementTranspose_C",NULL);CHKERRQ(ierr);
857*6dba178dSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorFactorizeSchurComplement_C",NULL);CHKERRQ(ierr);
858*6dba178dSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSchurComplementSetSym_C",NULL);CHKERRQ(ierr);
859bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetIcntl_C",NULL);CHKERRQ(ierr);
860bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetIcntl_C",NULL);CHKERRQ(ierr);
861bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetCntl_C",NULL);CHKERRQ(ierr);
862bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetCntl_C",NULL);CHKERRQ(ierr);
863ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfo_C",NULL);CHKERRQ(ierr);
864ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfog_C",NULL);CHKERRQ(ierr);
865ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfo_C",NULL);CHKERRQ(ierr);
866ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfog_C",NULL);CHKERRQ(ierr);
867397b6df1SKris Buschelman   PetscFunctionReturn(0);
868397b6df1SKris Buschelman }
869397b6df1SKris Buschelman 
870397b6df1SKris Buschelman #undef __FUNCT__
871f6c57405SHong Zhang #define __FUNCT__ "MatSolve_MUMPS"
872b24902e0SBarry Smith PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x)
873b24902e0SBarry Smith {
874a5e57a09SHong Zhang   Mat_MUMPS        *mumps=(Mat_MUMPS*)A->spptr;
875d54de34fSKris Buschelman   PetscScalar      *array;
87667877ebaSShri Abhyankar   Vec              b_seq;
877329ec9b3SHong Zhang   IS               is_iden,is_petsc;
878dfbe8321SBarry Smith   PetscErrorCode   ierr;
879329ec9b3SHong Zhang   PetscInt         i;
880cc86f929SStefano Zampini   PetscBool        second_solve = PETSC_FALSE;
881883f2eb9SBarry Smith   static PetscBool cite1 = PETSC_FALSE,cite2 = PETSC_FALSE;
882397b6df1SKris Buschelman 
883397b6df1SKris Buschelman   PetscFunctionBegin;
884883f2eb9SBarry 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);
885883f2eb9SBarry 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);
886a5e57a09SHong Zhang   mumps->id.nrhs = 1;
887a5e57a09SHong Zhang   b_seq          = mumps->b_seq;
888a5e57a09SHong Zhang   if (mumps->size > 1) {
889329ec9b3SHong Zhang     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
890a5e57a09SHong Zhang     ierr = VecScatterBegin(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
891a5e57a09SHong Zhang     ierr = VecScatterEnd(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
892a5e57a09SHong Zhang     if (!mumps->myid) {ierr = VecGetArray(b_seq,&array);CHKERRQ(ierr);}
893397b6df1SKris Buschelman   } else {  /* size == 1 */
894397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
895397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
896397b6df1SKris Buschelman   }
897a5e57a09SHong Zhang   if (!mumps->myid) { /* define rhs on the host */
898a5e57a09SHong Zhang     mumps->id.nrhs = 1;
899940cd9d6SSatish Balay     mumps->id.rhs = (MumpsScalar*)array;
900397b6df1SKris Buschelman   }
901397b6df1SKris Buschelman 
902cc86f929SStefano Zampini   /*
903cc86f929SStefano Zampini      handle condensation step of Schur complement (if any)
904cc86f929SStefano Zampini      We set by default ICNTL(26) == -1 when Schur indices have been provided by the user.
905cc86f929SStefano Zampini      According to MUMPS (5.0.0) manual, any value should be harmful during the factorization phase
906cc86f929SStefano Zampini      Unless the user provides a valid value for ICNTL(26), MatSolve and MatMatSolve routines solve the full system.
907cc86f929SStefano Zampini      This requires an extra call to PetscMUMPS_c and the computation of the factors for S
908cc86f929SStefano Zampini   */
909cc86f929SStefano Zampini   if (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2) {
910cc86f929SStefano Zampini     second_solve = PETSC_TRUE;
911b8f61ee1SStefano Zampini     ierr = MatMumpsHandleSchur_Private(mumps,PETSC_FALSE);CHKERRQ(ierr);
912cc86f929SStefano Zampini   }
913397b6df1SKris Buschelman   /* solve phase */
914329ec9b3SHong Zhang   /*-------------*/
915a5e57a09SHong Zhang   mumps->id.job = JOB_SOLVE;
916a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
917a5e57a09SHong 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));
918397b6df1SKris Buschelman 
919b5fa320bSStefano Zampini   /* handle expansion step of Schur complement (if any) */
920cc86f929SStefano Zampini   if (second_solve) {
921b8f61ee1SStefano Zampini     ierr = MatMumpsHandleSchur_Private(mumps,PETSC_TRUE);CHKERRQ(ierr);
922cc86f929SStefano Zampini   }
923b5fa320bSStefano Zampini 
924a5e57a09SHong Zhang   if (mumps->size > 1) { /* convert mumps distributed solution to petsc mpi x */
925a5e57a09SHong Zhang     if (mumps->scat_sol && mumps->ICNTL9_pre != mumps->id.ICNTL(9)) {
926a5e57a09SHong Zhang       /* when id.ICNTL(9) changes, the contents of lsol_loc may change (not its size, lsol_loc), recreates scat_sol */
927a5e57a09SHong Zhang       ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
928397b6df1SKris Buschelman     }
929a5e57a09SHong Zhang     if (!mumps->scat_sol) { /* create scatter scat_sol */
930a5e57a09SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,mumps->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */
931a5e57a09SHong Zhang       for (i=0; i<mumps->id.lsol_loc; i++) {
932a5e57a09SHong Zhang         mumps->id.isol_loc[i] -= 1; /* change Fortran style to C style */
933a5e57a09SHong Zhang       }
934a5e57a09SHong Zhang       ierr = ISCreateGeneral(PETSC_COMM_SELF,mumps->id.lsol_loc,mumps->id.isol_loc,PETSC_COPY_VALUES,&is_petsc);CHKERRQ(ierr);  /* to */
935a5e57a09SHong Zhang       ierr = VecScatterCreate(mumps->x_seq,is_iden,x,is_petsc,&mumps->scat_sol);CHKERRQ(ierr);
9366bf464f9SBarry Smith       ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
9376bf464f9SBarry Smith       ierr = ISDestroy(&is_petsc);CHKERRQ(ierr);
9382205254eSKarl Rupp 
939a5e57a09SHong Zhang       mumps->ICNTL9_pre = mumps->id.ICNTL(9); /* save current value of id.ICNTL(9) */
940397b6df1SKris Buschelman     }
941a5e57a09SHong Zhang 
942a5e57a09SHong Zhang     ierr = VecScatterBegin(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
943a5e57a09SHong Zhang     ierr = VecScatterEnd(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
944329ec9b3SHong Zhang   }
945397b6df1SKris Buschelman   PetscFunctionReturn(0);
946397b6df1SKris Buschelman }
947397b6df1SKris Buschelman 
94851d5961aSHong Zhang #undef __FUNCT__
94951d5961aSHong Zhang #define __FUNCT__ "MatSolveTranspose_MUMPS"
95051d5961aSHong Zhang PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x)
95151d5961aSHong Zhang {
952a5e57a09SHong Zhang   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
95351d5961aSHong Zhang   PetscErrorCode ierr;
95451d5961aSHong Zhang 
95551d5961aSHong Zhang   PetscFunctionBegin;
956a5e57a09SHong Zhang   mumps->id.ICNTL(9) = 0;
9570ad0caddSJed Brown   ierr = MatSolve_MUMPS(A,b,x);CHKERRQ(ierr);
958a5e57a09SHong Zhang   mumps->id.ICNTL(9) = 1;
95951d5961aSHong Zhang   PetscFunctionReturn(0);
96051d5961aSHong Zhang }
96151d5961aSHong Zhang 
962e0b74bf9SHong Zhang #undef __FUNCT__
963e0b74bf9SHong Zhang #define __FUNCT__ "MatMatSolve_MUMPS"
964e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X)
965e0b74bf9SHong Zhang {
966bda8bf91SBarry Smith   PetscErrorCode ierr;
967bda8bf91SBarry Smith   PetscBool      flg;
9684e34a73bSHong Zhang   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
969334c5f61SHong Zhang   PetscInt       i,nrhs,M;
9702cd7d884SHong Zhang   PetscScalar    *array,*bray;
971bda8bf91SBarry Smith 
972e0b74bf9SHong Zhang   PetscFunctionBegin;
9730298fd71SBarry Smith   ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
974801fbe65SHong Zhang   if (!flg) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
9750298fd71SBarry Smith   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
976801fbe65SHong Zhang   if (!flg) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
977801fbe65SHong 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");
9784e34a73bSHong Zhang 
9792cd7d884SHong Zhang   ierr = MatGetSize(B,&M,&nrhs);CHKERRQ(ierr);
980334c5f61SHong Zhang   mumps->id.nrhs = nrhs;
981334c5f61SHong Zhang   mumps->id.lrhs = M;
9824e34a73bSHong Zhang 
9832cd7d884SHong Zhang   if (mumps->size == 1) {
9842cd7d884SHong Zhang     /* copy B to X */
9852cd7d884SHong Zhang     ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr);
9862cd7d884SHong Zhang     ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr);
9876444a565SStefano Zampini     ierr = PetscMemcpy(array,bray,M*nrhs*sizeof(PetscScalar));CHKERRQ(ierr);
9882cd7d884SHong Zhang     ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr);
989940cd9d6SSatish Balay     mumps->id.rhs = (MumpsScalar*)array;
990b5fa320bSStefano Zampini     /* handle condensation step of Schur complement (if any) */
991b8f61ee1SStefano Zampini     ierr = MatMumpsHandleSchur_Private(mumps,PETSC_FALSE);CHKERRQ(ierr);
992801fbe65SHong Zhang 
9932cd7d884SHong Zhang     /* solve phase */
9942cd7d884SHong Zhang     /*-------------*/
9952cd7d884SHong Zhang     mumps->id.job = JOB_SOLVE;
9962cd7d884SHong Zhang     PetscMUMPS_c(&mumps->id);
9972cd7d884SHong 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));
998b5fa320bSStefano Zampini 
999b5fa320bSStefano Zampini     /* handle expansion step of Schur complement (if any) */
1000b8f61ee1SStefano Zampini     ierr = MatMumpsHandleSchur_Private(mumps,PETSC_TRUE);CHKERRQ(ierr);
10012cd7d884SHong Zhang     ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr);
1002334c5f61SHong Zhang   } else {  /*--------- parallel case --------*/
100371aed81dSHong Zhang     PetscInt       lsol_loc,nlsol_loc,*isol_loc,*idx,*iidx,*idxx,*isol_loc_save;
10041070efccSSatish Balay     MumpsScalar    *sol_loc,*sol_loc_save;
1005801fbe65SHong Zhang     IS             is_to,is_from;
1006334c5f61SHong Zhang     PetscInt       k,proc,j,m;
1007801fbe65SHong Zhang     const PetscInt *rstart;
1008334c5f61SHong Zhang     Vec            v_mpi,b_seq,x_seq;
1009334c5f61SHong Zhang     VecScatter     scat_rhs,scat_sol;
1010801fbe65SHong Zhang 
1011801fbe65SHong Zhang     /* create x_seq to hold local solution */
101271aed81dSHong Zhang     isol_loc_save = mumps->id.isol_loc; /* save it for MatSovle() */
101371aed81dSHong Zhang     sol_loc_save  = mumps->id.sol_loc;
1014801fbe65SHong Zhang 
101571aed81dSHong Zhang     lsol_loc  = mumps->id.INFO(23);
101671aed81dSHong Zhang     nlsol_loc = nrhs*lsol_loc;     /* length of sol_loc */
101771aed81dSHong Zhang     ierr = PetscMalloc2(nlsol_loc,&sol_loc,nlsol_loc,&isol_loc);CHKERRQ(ierr);
1018940cd9d6SSatish Balay     mumps->id.sol_loc = (MumpsScalar*)sol_loc;
1019801fbe65SHong Zhang     mumps->id.isol_loc = isol_loc;
1020801fbe65SHong Zhang 
10211070efccSSatish Balay     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,nlsol_loc,(PetscScalar*)sol_loc,&x_seq);CHKERRQ(ierr);
10222cd7d884SHong Zhang 
102374f0fcc7SHong Zhang     /* copy rhs matrix B into vector v_mpi */
1024334c5f61SHong Zhang     ierr = MatGetLocalSize(B,&m,NULL);CHKERRQ(ierr);
1025801fbe65SHong Zhang     ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr);
102674f0fcc7SHong Zhang     ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)B),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi);CHKERRQ(ierr);
1027801fbe65SHong Zhang     ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr);
1028801fbe65SHong Zhang 
1029334c5f61SHong Zhang     /* scatter v_mpi to b_seq because MUMPS only supports centralized rhs */
103074f0fcc7SHong Zhang     /* idx: maps from k-th index of v_mpi to (i,j)-th global entry of B;
1031801fbe65SHong Zhang       iidx: inverse of idx, will be used by scattering xx_seq -> X       */
1032801fbe65SHong Zhang     ierr = PetscMalloc2(nrhs*M,&idx,nrhs*M,&iidx);CHKERRQ(ierr);
1033801fbe65SHong Zhang     ierr = MatGetOwnershipRanges(B,&rstart);CHKERRQ(ierr);
1034801fbe65SHong Zhang     k = 0;
1035801fbe65SHong Zhang     for (proc=0; proc<mumps->size; proc++){
1036801fbe65SHong Zhang       for (j=0; j<nrhs; j++){
1037801fbe65SHong Zhang         for (i=rstart[proc]; i<rstart[proc+1]; i++){
1038801fbe65SHong Zhang           iidx[j*M + i] = k;
1039801fbe65SHong Zhang           idx[k++]      = j*M + i;
1040801fbe65SHong Zhang         }
1041801fbe65SHong Zhang       }
10422cd7d884SHong Zhang     }
10432cd7d884SHong Zhang 
1044801fbe65SHong Zhang     if (!mumps->myid) {
1045334c5f61SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,nrhs*M,&b_seq);CHKERRQ(ierr);
1046801fbe65SHong Zhang       ierr = ISCreateGeneral(PETSC_COMM_SELF,nrhs*M,idx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr);
1047801fbe65SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,nrhs*M,0,1,&is_from);CHKERRQ(ierr);
1048801fbe65SHong Zhang     } else {
1049334c5f61SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&b_seq);CHKERRQ(ierr);
1050801fbe65SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_to);CHKERRQ(ierr);
1051801fbe65SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_from);CHKERRQ(ierr);
1052801fbe65SHong Zhang     }
1053334c5f61SHong Zhang     ierr = VecScatterCreate(v_mpi,is_from,b_seq,is_to,&scat_rhs);CHKERRQ(ierr);
1054334c5f61SHong Zhang     ierr = VecScatterBegin(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1055801fbe65SHong Zhang     ierr = ISDestroy(&is_to);CHKERRQ(ierr);
1056801fbe65SHong Zhang     ierr = ISDestroy(&is_from);CHKERRQ(ierr);
1057334c5f61SHong Zhang     ierr = VecScatterEnd(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1058801fbe65SHong Zhang 
1059801fbe65SHong Zhang     if (!mumps->myid) { /* define rhs on the host */
1060334c5f61SHong Zhang       ierr = VecGetArray(b_seq,&bray);CHKERRQ(ierr);
1061940cd9d6SSatish Balay       mumps->id.rhs = (MumpsScalar*)bray;
1062334c5f61SHong Zhang       ierr = VecRestoreArray(b_seq,&bray);CHKERRQ(ierr);
1063801fbe65SHong Zhang     }
1064801fbe65SHong Zhang 
1065801fbe65SHong Zhang     /* solve phase */
1066801fbe65SHong Zhang     /*-------------*/
1067801fbe65SHong Zhang     mumps->id.job = JOB_SOLVE;
1068801fbe65SHong Zhang     PetscMUMPS_c(&mumps->id);
1069801fbe65SHong 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));
1070801fbe65SHong Zhang 
1071334c5f61SHong Zhang     /* scatter mumps distributed solution to petsc vector v_mpi, which shares local arrays with solution matrix X */
107274f0fcc7SHong Zhang     ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr);
107374f0fcc7SHong Zhang     ierr = VecPlaceArray(v_mpi,array);CHKERRQ(ierr);
1074801fbe65SHong Zhang 
1075334c5f61SHong Zhang     /* create scatter scat_sol */
107671aed81dSHong Zhang     ierr = PetscMalloc1(nlsol_loc,&idxx);CHKERRQ(ierr);
107771aed81dSHong Zhang     ierr = ISCreateStride(PETSC_COMM_SELF,nlsol_loc,0,1,&is_from);CHKERRQ(ierr);
107871aed81dSHong Zhang     for (i=0; i<lsol_loc; i++) {
1079334c5f61SHong Zhang       isol_loc[i] -= 1; /* change Fortran style to C style */
1080334c5f61SHong Zhang       idxx[i] = iidx[isol_loc[i]];
1081801fbe65SHong Zhang       for (j=1; j<nrhs; j++){
1082334c5f61SHong Zhang         idxx[j*lsol_loc+i] = iidx[isol_loc[i]+j*M];
1083801fbe65SHong Zhang       }
1084801fbe65SHong Zhang     }
108571aed81dSHong Zhang     ierr = ISCreateGeneral(PETSC_COMM_SELF,nlsol_loc,idxx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr);
1086334c5f61SHong Zhang     ierr = VecScatterCreate(x_seq,is_from,v_mpi,is_to,&scat_sol);CHKERRQ(ierr);
1087334c5f61SHong Zhang     ierr = VecScatterBegin(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1088801fbe65SHong Zhang     ierr = ISDestroy(&is_from);CHKERRQ(ierr);
1089801fbe65SHong Zhang     ierr = ISDestroy(&is_to);CHKERRQ(ierr);
1090334c5f61SHong Zhang     ierr = VecScatterEnd(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1091801fbe65SHong Zhang     ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr);
109271aed81dSHong Zhang 
109371aed81dSHong Zhang     /* free spaces */
109471aed81dSHong Zhang     mumps->id.sol_loc = sol_loc_save;
109571aed81dSHong Zhang     mumps->id.isol_loc = isol_loc_save;
109671aed81dSHong Zhang 
109771aed81dSHong Zhang     ierr = PetscFree2(sol_loc,isol_loc);CHKERRQ(ierr);
1098801fbe65SHong Zhang     ierr = PetscFree2(idx,iidx);CHKERRQ(ierr);
1099801fbe65SHong Zhang     ierr = PetscFree(idxx);CHKERRQ(ierr);
110071aed81dSHong Zhang     ierr = VecDestroy(&x_seq);CHKERRQ(ierr);
110174f0fcc7SHong Zhang     ierr = VecDestroy(&v_mpi);CHKERRQ(ierr);
1102334c5f61SHong Zhang     ierr = VecDestroy(&b_seq);CHKERRQ(ierr);
1103334c5f61SHong Zhang     ierr = VecScatterDestroy(&scat_rhs);CHKERRQ(ierr);
1104334c5f61SHong Zhang     ierr = VecScatterDestroy(&scat_sol);CHKERRQ(ierr);
1105801fbe65SHong Zhang   }
1106e0b74bf9SHong Zhang   PetscFunctionReturn(0);
1107e0b74bf9SHong Zhang }
1108e0b74bf9SHong Zhang 
1109ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX)
1110a58c3f20SHong Zhang /*
1111a58c3f20SHong Zhang   input:
1112a58c3f20SHong Zhang    F:        numeric factor
1113a58c3f20SHong Zhang   output:
1114a58c3f20SHong Zhang    nneg:     total number of negative pivots
1115a58c3f20SHong Zhang    nzero:    0
1116a58c3f20SHong Zhang    npos:     (global dimension of F) - nneg
1117a58c3f20SHong Zhang */
1118a58c3f20SHong Zhang 
1119a58c3f20SHong Zhang #undef __FUNCT__
1120a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
1121dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
1122a58c3f20SHong Zhang {
1123a5e57a09SHong Zhang   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
1124dfbe8321SBarry Smith   PetscErrorCode ierr;
1125c1490034SHong Zhang   PetscMPIInt    size;
1126a58c3f20SHong Zhang 
1127a58c3f20SHong Zhang   PetscFunctionBegin;
1128ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)F),&size);CHKERRQ(ierr);
1129bcb30aebSHong 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 */
1130a5e57a09SHong 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));
1131ed85ac9fSHong Zhang 
1132710ac8efSHong Zhang   if (nneg) *nneg = mumps->id.INFOG(12);
1133ed85ac9fSHong Zhang   if (nzero || npos) {
1134ed85ac9fSHong 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");
1135710ac8efSHong Zhang     if (nzero) *nzero = mumps->id.INFOG(28);
1136710ac8efSHong Zhang     if (npos) *npos   = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28));
1137a58c3f20SHong Zhang   }
1138a58c3f20SHong Zhang   PetscFunctionReturn(0);
1139a58c3f20SHong Zhang }
1140ace3df97SHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */
1141a58c3f20SHong Zhang 
1142397b6df1SKris Buschelman #undef __FUNCT__
1143f6c57405SHong Zhang #define __FUNCT__ "MatFactorNumeric_MUMPS"
11440481f469SBarry Smith PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
1145af281ebdSHong Zhang {
1146a5e57a09SHong Zhang   Mat_MUMPS      *mumps =(Mat_MUMPS*)(F)->spptr;
11476849ba73SBarry Smith   PetscErrorCode ierr;
1148e09efc27SHong Zhang   Mat            F_diag;
1149ace3abfcSBarry Smith   PetscBool      isMPIAIJ;
1150397b6df1SKris Buschelman 
1151397b6df1SKris Buschelman   PetscFunctionBegin;
1152a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1153397b6df1SKris Buschelman 
1154397b6df1SKris Buschelman   /* numerical factorization phase */
1155329ec9b3SHong Zhang   /*-------------------------------*/
1156a5e57a09SHong Zhang   mumps->id.job = JOB_FACTNUMERIC;
11574e34a73bSHong Zhang   if (!mumps->id.ICNTL(18)) { /* A is centralized */
1158a5e57a09SHong Zhang     if (!mumps->myid) {
1159940cd9d6SSatish Balay       mumps->id.a = (MumpsScalar*)mumps->val;
1160397b6df1SKris Buschelman     }
1161397b6df1SKris Buschelman   } else {
1162940cd9d6SSatish Balay     mumps->id.a_loc = (MumpsScalar*)mumps->val;
1163397b6df1SKris Buschelman   }
1164a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
1165a5e57a09SHong Zhang   if (mumps->id.INFOG(1) < 0) {
1166151787a6SHong Zhang     if (mumps->id.INFO(1) == -13) {
1167151787a6SHong Zhang       if (mumps->id.INFO(2) < 0) {
1168151787a6SHong Zhang         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",-mumps->id.INFO(2));
1169151787a6SHong Zhang       } else {
1170151787a6SHong Zhang         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d bytes\n",mumps->id.INFO(2));
1171151787a6SHong Zhang       }
1172151787a6SHong Zhang     } else SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: INFO(1)=%d, INFO(2)=%d\n",mumps->id.INFO(1),mumps->id.INFO(2));
1173397b6df1SKris Buschelman   }
1174a5e57a09SHong 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));
1175397b6df1SKris Buschelman 
1176dcd589f8SShri Abhyankar   (F)->assembled        = PETSC_TRUE;
1177a5e57a09SHong Zhang   mumps->matstruc       = SAME_NONZERO_PATTERN;
1178b5fa320bSStefano Zampini   mumps->schur_factored = PETSC_FALSE;
117959ac8732SStefano Zampini   mumps->schur_inverted = PETSC_FALSE;
118067877ebaSShri Abhyankar 
1181066565c5SStefano Zampini   /* just to be sure that ICNTL(19) value returned by a call from MatMumpsGetIcntl is always consistent */
1182066565c5SStefano Zampini   if (!mumps->sym && mumps->id.ICNTL(19) && mumps->id.ICNTL(19) != 1) mumps->id.ICNTL(19) = 3;
1183066565c5SStefano Zampini 
1184a5e57a09SHong Zhang   if (mumps->size > 1) {
118567877ebaSShri Abhyankar     PetscInt    lsol_loc;
118667877ebaSShri Abhyankar     PetscScalar *sol_loc;
11872205254eSKarl Rupp 
1188c2093ab7SHong Zhang     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr);
1189c2093ab7SHong Zhang     if (isMPIAIJ) F_diag = ((Mat_MPIAIJ*)(F)->data)->A;
1190c2093ab7SHong Zhang     else F_diag = ((Mat_MPISBAIJ*)(F)->data)->A;
1191c2093ab7SHong Zhang     F_diag->assembled = PETSC_TRUE;
1192c2093ab7SHong Zhang 
1193c2093ab7SHong Zhang     /* distributed solution; Create x_seq=sol_loc for repeated use */
1194c2093ab7SHong Zhang     if (mumps->x_seq) {
1195c2093ab7SHong Zhang       ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
1196c2093ab7SHong Zhang       ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr);
1197c2093ab7SHong Zhang       ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr);
1198c2093ab7SHong Zhang     }
1199a5e57a09SHong Zhang     lsol_loc = mumps->id.INFO(23); /* length of sol_loc */
1200dcca6d9dSJed Brown     ierr = PetscMalloc2(lsol_loc,&sol_loc,lsol_loc,&mumps->id.isol_loc);CHKERRQ(ierr);
1201a5e57a09SHong Zhang     mumps->id.lsol_loc = lsol_loc;
1202940cd9d6SSatish Balay     mumps->id.sol_loc = (MumpsScalar*)sol_loc;
1203a5e57a09SHong Zhang     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq);CHKERRQ(ierr);
120467877ebaSShri Abhyankar   }
1205397b6df1SKris Buschelman   PetscFunctionReturn(0);
1206397b6df1SKris Buschelman }
1207397b6df1SKris Buschelman 
12089a2535b5SHong Zhang /* Sets MUMPS options from the options database */
1209dcd589f8SShri Abhyankar #undef __FUNCT__
12109a2535b5SHong Zhang #define __FUNCT__ "PetscSetMUMPSFromOptions"
12119a2535b5SHong Zhang PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A)
1212dcd589f8SShri Abhyankar {
12139a2535b5SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
1214dcd589f8SShri Abhyankar   PetscErrorCode ierr;
1215b34f08ffSHong Zhang   PetscInt       icntl,info[40],i,ninfo=40;
1216ace3abfcSBarry Smith   PetscBool      flg;
1217dcd589f8SShri Abhyankar 
1218dcd589f8SShri Abhyankar   PetscFunctionBegin;
1219ce94432eSBarry Smith   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
12209a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);CHKERRQ(ierr);
12219a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(1) = icntl;
12229a2535b5SHong 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);
12239a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(2) = icntl;
12249a2535b5SHong 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);
12259a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(3) = icntl;
1226dcd589f8SShri Abhyankar 
12279a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",mumps->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
12289a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(4) = icntl;
12299a2535b5SHong Zhang   if (mumps->id.ICNTL(4) || PetscLogPrintInfo) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */
12309a2535b5SHong Zhang 
1231d341cd04SHong 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);
12329a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(6) = icntl;
12339a2535b5SHong Zhang 
1234d341cd04SHong 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);
1235dcd589f8SShri Abhyankar   if (flg) {
12362205254eSKarl 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");
12372205254eSKarl Rupp     else mumps->id.ICNTL(7) = icntl;
1238dcd589f8SShri Abhyankar   }
1239e0b74bf9SHong Zhang 
12400298fd71SBarry 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);
1241d341cd04SHong 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() */
12420298fd71SBarry 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);
1243d341cd04SHong 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);
1244d341cd04SHong 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);
1245d341cd04SHong 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);
1246d341cd04SHong 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);
1247d341cd04SHong 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);
124859ac8732SStefano Zampini   if (mumps->id.ICNTL(19) <= 0 || mumps->id.ICNTL(19) > 3) { /* reset any schur data (if any) */
124959ac8732SStefano Zampini     ierr = MatMumpsResetSchur_Private(mumps);CHKERRQ(ierr);
125059ac8732SStefano Zampini   }
12514e34a73bSHong 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 */
1252d341cd04SHong 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 */
12539a2535b5SHong Zhang 
1254d341cd04SHong 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);
12550298fd71SBarry 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);
12560298fd71SBarry 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);
12579a2535b5SHong Zhang   if (mumps->id.ICNTL(24)) {
12589a2535b5SHong Zhang     mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */
1259d7ebd59bSHong Zhang   }
1260d7ebd59bSHong Zhang 
1261d341cd04SHong 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);
1262d341cd04SHong 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);
12632cd7d884SHong 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);
12640298fd71SBarry 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);
1265d341cd04SHong 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);
12660298fd71SBarry 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);
1267d341cd04SHong 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);
12684e34a73bSHong 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 */
12690298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),NULL);CHKERRQ(ierr);
1270dcd589f8SShri Abhyankar 
12710298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),NULL);CHKERRQ(ierr);
12720298fd71SBarry 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);
12730298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),NULL);CHKERRQ(ierr);
12740298fd71SBarry 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);
12750298fd71SBarry 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);
1276e5bb22a1SHong Zhang 
12770298fd71SBarry Smith   ierr = PetscOptionsString("-mat_mumps_ooc_tmpdir", "out of core directory", "None", mumps->id.ooc_tmpdir, mumps->id.ooc_tmpdir, 256, NULL);
1278b34f08ffSHong Zhang 
127916d797efSHong Zhang   ierr = PetscOptionsIntArray("-mat_mumps_view_info","request INFO local to each processor","",info,&ninfo,NULL);CHKERRQ(ierr);
1280b34f08ffSHong Zhang   if (ninfo) {
1281b34f08ffSHong Zhang     if (ninfo > 40) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"number of INFO %d must <= 40\n",ninfo);
1282b34f08ffSHong Zhang     ierr = PetscMalloc1(ninfo,&mumps->info);CHKERRQ(ierr);
1283b34f08ffSHong Zhang     mumps->ninfo = ninfo;
1284b34f08ffSHong Zhang     for (i=0; i<ninfo; i++) {
1285b34f08ffSHong Zhang       if (info[i] < 0 || info[i]>40) {
1286b34f08ffSHong Zhang         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index of INFO %d must between 1 and 40\n",ninfo);
1287b34f08ffSHong Zhang       } else {
1288b34f08ffSHong Zhang         mumps->info[i] = info[i];
1289b34f08ffSHong Zhang       }
1290b34f08ffSHong Zhang     }
1291b34f08ffSHong Zhang   }
1292b34f08ffSHong Zhang 
1293dcd589f8SShri Abhyankar   PetscOptionsEnd();
1294dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
1295dcd589f8SShri Abhyankar }
1296dcd589f8SShri Abhyankar 
1297dcd589f8SShri Abhyankar #undef __FUNCT__
1298dcd589f8SShri Abhyankar #define __FUNCT__ "PetscInitializeMUMPS"
1299f697e70eSHong Zhang PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps)
1300dcd589f8SShri Abhyankar {
1301dcd589f8SShri Abhyankar   PetscErrorCode ierr;
1302dcd589f8SShri Abhyankar 
1303dcd589f8SShri Abhyankar   PetscFunctionBegin;
1304ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)A), &mumps->myid);
1305ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&mumps->size);CHKERRQ(ierr);
1306ce94432eSBarry Smith   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(mumps->comm_mumps));CHKERRQ(ierr);
13072205254eSKarl Rupp 
1308f697e70eSHong Zhang   mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps);
1309f697e70eSHong Zhang 
1310f697e70eSHong Zhang   mumps->id.job = JOB_INIT;
1311f697e70eSHong Zhang   mumps->id.par = 1;  /* host participates factorizaton and solve */
1312f697e70eSHong Zhang   mumps->id.sym = mumps->sym;
13132907cef9SHong Zhang   PetscMUMPS_c(&mumps->id);
1314f697e70eSHong Zhang 
13150298fd71SBarry Smith   mumps->scat_rhs     = NULL;
13160298fd71SBarry Smith   mumps->scat_sol     = NULL;
13179a2535b5SHong Zhang 
131870544d5fSHong Zhang   /* set PETSc-MUMPS default options - override MUMPS default */
13199a2535b5SHong Zhang   mumps->id.ICNTL(3) = 0;
13209a2535b5SHong Zhang   mumps->id.ICNTL(4) = 0;
13219a2535b5SHong Zhang   if (mumps->size == 1) {
13229a2535b5SHong Zhang     mumps->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
13239a2535b5SHong Zhang   } else {
13249a2535b5SHong Zhang     mumps->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
13254e34a73bSHong Zhang     mumps->id.ICNTL(20) = 0;   /* rhs is in dense format */
132670544d5fSHong Zhang     mumps->id.ICNTL(21) = 1;   /* distributed solution */
13279a2535b5SHong Zhang   }
13286444a565SStefano Zampini 
13296444a565SStefano Zampini   /* schur */
13306444a565SStefano Zampini   mumps->id.size_schur      = 0;
13316444a565SStefano Zampini   mumps->id.listvar_schur   = NULL;
13326444a565SStefano Zampini   mumps->id.schur           = NULL;
1333b5fa320bSStefano Zampini   mumps->sizeredrhs         = 0;
1334b5fa320bSStefano Zampini   mumps->schur_pivots       = NULL;
1335b5fa320bSStefano Zampini   mumps->schur_work         = NULL;
133659ac8732SStefano Zampini   mumps->schur_sol          = NULL;
133759ac8732SStefano Zampini   mumps->schur_sizesol      = 0;
133859ac8732SStefano Zampini   mumps->schur_factored     = PETSC_FALSE;
133959ac8732SStefano Zampini   mumps->schur_inverted     = PETSC_FALSE;
1340a0b0af32SStefano Zampini   mumps->schur_sym          = mumps->id.sym;
1341dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
1342dcd589f8SShri Abhyankar }
1343dcd589f8SShri Abhyankar 
1344a5e57a09SHong 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 */
1345397b6df1SKris Buschelman #undef __FUNCT__
1346f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
13470481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1348b24902e0SBarry Smith {
1349a5e57a09SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
1350dcd589f8SShri Abhyankar   PetscErrorCode ierr;
135167877ebaSShri Abhyankar   Vec            b;
135267877ebaSShri Abhyankar   IS             is_iden;
135367877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1354397b6df1SKris Buschelman 
1355397b6df1SKris Buschelman   PetscFunctionBegin;
1356a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1357dcd589f8SShri Abhyankar 
13589a2535b5SHong Zhang   /* Set MUMPS options from the options database */
13599a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1360dcd589f8SShri Abhyankar 
1361a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1362dcd589f8SShri Abhyankar 
136367877ebaSShri Abhyankar   /* analysis phase */
136467877ebaSShri Abhyankar   /*----------------*/
1365a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1366a5e57a09SHong Zhang   mumps->id.n   = M;
1367a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
136867877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1369a5e57a09SHong Zhang     if (!mumps->myid) {
1370a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1371a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
1372940cd9d6SSatish Balay         mumps->id.a = (MumpsScalar*)mumps->val;
137367877ebaSShri Abhyankar       }
1374a5e57a09SHong Zhang       if (mumps->id.ICNTL(7) == 1) { /* use user-provide matrix ordering - assuming r = c ordering */
13755248a706SHong Zhang         /*
13765248a706SHong Zhang         PetscBool      flag;
13775248a706SHong Zhang         ierr = ISEqual(r,c,&flag);CHKERRQ(ierr);
13785248a706SHong Zhang         if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"row_perm != col_perm");
13795248a706SHong Zhang         ierr = ISView(r,PETSC_VIEWER_STDOUT_SELF);
13805248a706SHong Zhang          */
1381a5e57a09SHong Zhang         if (!mumps->myid) {
1382e0b74bf9SHong Zhang           const PetscInt *idx;
1383e0b74bf9SHong Zhang           PetscInt       i,*perm_in;
13842205254eSKarl Rupp 
1385785e854fSJed Brown           ierr = PetscMalloc1(M,&perm_in);CHKERRQ(ierr);
1386e0b74bf9SHong Zhang           ierr = ISGetIndices(r,&idx);CHKERRQ(ierr);
13872205254eSKarl Rupp 
1388a5e57a09SHong Zhang           mumps->id.perm_in = perm_in;
1389e0b74bf9SHong Zhang           for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */
1390e0b74bf9SHong Zhang           ierr = ISRestoreIndices(r,&idx);CHKERRQ(ierr);
1391e0b74bf9SHong Zhang         }
1392e0b74bf9SHong Zhang       }
139367877ebaSShri Abhyankar     }
139467877ebaSShri Abhyankar     break;
139567877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1396a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1397a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1398a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
1399940cd9d6SSatish Balay       mumps->id.a_loc = (MumpsScalar*)mumps->val;
140067877ebaSShri Abhyankar     }
140167877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1402a5e57a09SHong Zhang     if (!mumps->myid) {
14032cd7d884SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->rmap->N,&mumps->b_seq);CHKERRQ(ierr);
14042cd7d884SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->N,0,1,&is_iden);CHKERRQ(ierr);
140567877ebaSShri Abhyankar     } else {
1406a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
140767877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
140867877ebaSShri Abhyankar     }
14092a7a6963SBarry Smith     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1410a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
14116bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
14126bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
141367877ebaSShri Abhyankar     break;
141467877ebaSShri Abhyankar   }
1415a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
1416a5e57a09SHong Zhang   if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
141767877ebaSShri Abhyankar 
1418719d5645SBarry Smith   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1419dcd589f8SShri Abhyankar   F->ops->solve           = MatSolve_MUMPS;
142051d5961aSHong Zhang   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
14214e34a73bSHong Zhang   F->ops->matsolve        = MatMatSolve_MUMPS;
1422b24902e0SBarry Smith   PetscFunctionReturn(0);
1423b24902e0SBarry Smith }
1424b24902e0SBarry Smith 
1425450b117fSShri Abhyankar /* Note the Petsc r and c permutations are ignored */
1426450b117fSShri Abhyankar #undef __FUNCT__
1427450b117fSShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_BAIJMUMPS"
1428450b117fSShri Abhyankar PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1429450b117fSShri Abhyankar {
1430a5e57a09SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
1431dcd589f8SShri Abhyankar   PetscErrorCode ierr;
143267877ebaSShri Abhyankar   Vec            b;
143367877ebaSShri Abhyankar   IS             is_iden;
143467877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1435450b117fSShri Abhyankar 
1436450b117fSShri Abhyankar   PetscFunctionBegin;
1437a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1438dcd589f8SShri Abhyankar 
14399a2535b5SHong Zhang   /* Set MUMPS options from the options database */
14409a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1441dcd589f8SShri Abhyankar 
1442a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
144367877ebaSShri Abhyankar 
144467877ebaSShri Abhyankar   /* analysis phase */
144567877ebaSShri Abhyankar   /*----------------*/
1446a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1447a5e57a09SHong Zhang   mumps->id.n   = M;
1448a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
144967877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1450a5e57a09SHong Zhang     if (!mumps->myid) {
1451a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1452a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
1453940cd9d6SSatish Balay         mumps->id.a = (MumpsScalar*)mumps->val;
145467877ebaSShri Abhyankar       }
145567877ebaSShri Abhyankar     }
145667877ebaSShri Abhyankar     break;
145767877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1458a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1459a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1460a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
1461940cd9d6SSatish Balay       mumps->id.a_loc = (MumpsScalar*)mumps->val;
146267877ebaSShri Abhyankar     }
146367877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1464a5e57a09SHong Zhang     if (!mumps->myid) {
1465a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
146667877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
146767877ebaSShri Abhyankar     } else {
1468a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
146967877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
147067877ebaSShri Abhyankar     }
14712a7a6963SBarry Smith     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1472a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
14736bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
14746bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
147567877ebaSShri Abhyankar     break;
147667877ebaSShri Abhyankar   }
1477a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
1478a5e57a09SHong Zhang   if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
147967877ebaSShri Abhyankar 
1480450b117fSShri Abhyankar   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1481dcd589f8SShri Abhyankar   F->ops->solve           = MatSolve_MUMPS;
148251d5961aSHong Zhang   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1483450b117fSShri Abhyankar   PetscFunctionReturn(0);
1484450b117fSShri Abhyankar }
1485b24902e0SBarry Smith 
1486141f4205SHong Zhang /* Note the Petsc r permutation and factor info are ignored */
1487397b6df1SKris Buschelman #undef __FUNCT__
148867877ebaSShri Abhyankar #define __FUNCT__ "MatCholeskyFactorSymbolic_MUMPS"
148967877ebaSShri Abhyankar PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
1490b24902e0SBarry Smith {
1491a5e57a09SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
1492dcd589f8SShri Abhyankar   PetscErrorCode ierr;
149367877ebaSShri Abhyankar   Vec            b;
149467877ebaSShri Abhyankar   IS             is_iden;
149567877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1496397b6df1SKris Buschelman 
1497397b6df1SKris Buschelman   PetscFunctionBegin;
1498a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1499dcd589f8SShri Abhyankar 
15009a2535b5SHong Zhang   /* Set MUMPS options from the options database */
15019a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1502dcd589f8SShri Abhyankar 
1503a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1504dcd589f8SShri Abhyankar 
150567877ebaSShri Abhyankar   /* analysis phase */
150667877ebaSShri Abhyankar   /*----------------*/
1507a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1508a5e57a09SHong Zhang   mumps->id.n   = M;
1509a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
151067877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1511a5e57a09SHong Zhang     if (!mumps->myid) {
1512a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1513a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
1514940cd9d6SSatish Balay         mumps->id.a = (MumpsScalar*)mumps->val;
151567877ebaSShri Abhyankar       }
151667877ebaSShri Abhyankar     }
151767877ebaSShri Abhyankar     break;
151867877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1519a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1520a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1521a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
1522940cd9d6SSatish Balay       mumps->id.a_loc = (MumpsScalar*)mumps->val;
152367877ebaSShri Abhyankar     }
152467877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1525a5e57a09SHong Zhang     if (!mumps->myid) {
1526a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
152767877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
152867877ebaSShri Abhyankar     } else {
1529a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
153067877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
153167877ebaSShri Abhyankar     }
15322a7a6963SBarry Smith     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1533a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
15346bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
15356bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
153667877ebaSShri Abhyankar     break;
153767877ebaSShri Abhyankar   }
1538a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
1539a5e57a09SHong Zhang   if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
154067877ebaSShri Abhyankar 
15412792810eSHong Zhang   F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
1542dcd589f8SShri Abhyankar   F->ops->solve                 = MatSolve_MUMPS;
154351d5961aSHong Zhang   F->ops->solvetranspose        = MatSolve_MUMPS;
15444e34a73bSHong Zhang   F->ops->matsolve              = MatMatSolve_MUMPS;
15454e34a73bSHong Zhang #if defined(PETSC_USE_COMPLEX)
15460298fd71SBarry Smith   F->ops->getinertia = NULL;
15474e34a73bSHong Zhang #else
15484e34a73bSHong Zhang   F->ops->getinertia = MatGetInertia_SBAIJMUMPS;
1549db4efbfdSBarry Smith #endif
1550b24902e0SBarry Smith   PetscFunctionReturn(0);
1551b24902e0SBarry Smith }
1552b24902e0SBarry Smith 
1553397b6df1SKris Buschelman #undef __FUNCT__
155464e6c443SBarry Smith #define __FUNCT__ "MatView_MUMPS"
155564e6c443SBarry Smith PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
155674ed9c26SBarry Smith {
1557f6c57405SHong Zhang   PetscErrorCode    ierr;
155864e6c443SBarry Smith   PetscBool         iascii;
155964e6c443SBarry Smith   PetscViewerFormat format;
1560a5e57a09SHong Zhang   Mat_MUMPS         *mumps=(Mat_MUMPS*)A->spptr;
1561f6c57405SHong Zhang 
1562f6c57405SHong Zhang   PetscFunctionBegin;
156364e6c443SBarry Smith   /* check if matrix is mumps type */
156464e6c443SBarry Smith   if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0);
156564e6c443SBarry Smith 
1566251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
156764e6c443SBarry Smith   if (iascii) {
156864e6c443SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
156964e6c443SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO) {
157064e6c443SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
1571a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                   %d \n",mumps->id.sym);CHKERRQ(ierr);
1572a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):            %d \n",mumps->id.par);CHKERRQ(ierr);
1573a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):         %d \n",mumps->id.ICNTL(1));CHKERRQ(ierr);
1574a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg): %d \n",mumps->id.ICNTL(2));CHKERRQ(ierr);
1575a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):   %d \n",mumps->id.ICNTL(3));CHKERRQ(ierr);
1576a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):        %d \n",mumps->id.ICNTL(4));CHKERRQ(ierr);
1577a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):         %d \n",mumps->id.ICNTL(5));CHKERRQ(ierr);
1578a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):        %d \n",mumps->id.ICNTL(6));CHKERRQ(ierr);
1579a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (sequentia matrix ordering):%d \n",mumps->id.ICNTL(7));CHKERRQ(ierr);
1580a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scalling strategy):        %d \n",mumps->id.ICNTL(8));CHKERRQ(ierr);
1581a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements):  %d \n",mumps->id.ICNTL(10));CHKERRQ(ierr);
1582a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):          %d \n",mumps->id.ICNTL(11));CHKERRQ(ierr);
1583a5e57a09SHong Zhang       if (mumps->id.ICNTL(11)>0) {
1584a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(4) (inf norm of input mat):        %g\n",mumps->id.RINFOG(4));CHKERRQ(ierr);
1585a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(5) (inf norm of solution):         %g\n",mumps->id.RINFOG(5));CHKERRQ(ierr);
1586a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(6) (inf norm of residual):         %g\n",mumps->id.RINFOG(6));CHKERRQ(ierr);
1587a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8));CHKERRQ(ierr);
1588a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(9) (error estimate):               %g \n",mumps->id.RINFOG(9));CHKERRQ(ierr);
1589a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11));CHKERRQ(ierr);
1590f6c57405SHong Zhang       }
1591a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",mumps->id.ICNTL(12));CHKERRQ(ierr);
1592a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",mumps->id.ICNTL(13));CHKERRQ(ierr);
1593a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",mumps->id.ICNTL(14));CHKERRQ(ierr);
1594f6c57405SHong Zhang       /* ICNTL(15-17) not used */
1595a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",mumps->id.ICNTL(18));CHKERRQ(ierr);
1596a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Shur complement info):                       %d \n",mumps->id.ICNTL(19));CHKERRQ(ierr);
1597a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",mumps->id.ICNTL(20));CHKERRQ(ierr);
1598ca3dc52bSPierre Jolivet       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d \n",mumps->id.ICNTL(21));CHKERRQ(ierr);
1599a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",mumps->id.ICNTL(22));CHKERRQ(ierr);
1600a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",mumps->id.ICNTL(23));CHKERRQ(ierr);
1601c0165424SHong Zhang 
1602a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",mumps->id.ICNTL(24));CHKERRQ(ierr);
1603a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",mumps->id.ICNTL(25));CHKERRQ(ierr);
1604a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",mumps->id.ICNTL(26));CHKERRQ(ierr);
1605a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",mumps->id.ICNTL(27));CHKERRQ(ierr);
1606a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(28) (use parallel or sequential ordering):        %d \n",mumps->id.ICNTL(28));CHKERRQ(ierr);
1607a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(29) (parallel ordering):                          %d \n",mumps->id.ICNTL(29));CHKERRQ(ierr);
160842179a6aSHong Zhang 
1609a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(30) (user-specified set of entries in inv(A)):    %d \n",mumps->id.ICNTL(30));CHKERRQ(ierr);
1610a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(31) (factors is discarded in the solve phase):    %d \n",mumps->id.ICNTL(31));CHKERRQ(ierr);
1611a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(33) (compute determinant):                        %d \n",mumps->id.ICNTL(33));CHKERRQ(ierr);
1612f6c57405SHong Zhang 
1613a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",mumps->id.CNTL(1));CHKERRQ(ierr);
1614a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2));CHKERRQ(ierr);
1615ca3dc52bSPierre Jolivet       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",mumps->id.CNTL(3));CHKERRQ(ierr);
1616ca3dc52bSPierre Jolivet       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",mumps->id.CNTL(4));CHKERRQ(ierr);
1617a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",mumps->id.CNTL(5));CHKERRQ(ierr);
1618f6c57405SHong Zhang 
1619f6c57405SHong Zhang       /* infomation local to each processor */
162034ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);
16211575c14dSBarry Smith       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
1622a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %g \n",mumps->myid,mumps->id.RINFO(1));CHKERRQ(ierr);
162334ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
162434ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);
1625a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(2));CHKERRQ(ierr);
162634ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
162734ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);
1628a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(3));CHKERRQ(ierr);
162934ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1630f6c57405SHong Zhang 
163134ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr);
1632a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] %d \n",mumps->myid,mumps->id.INFO(15));CHKERRQ(ierr);
163334ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1634f6c57405SHong Zhang 
163534ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr);
1636a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(16));CHKERRQ(ierr);
163734ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1638f6c57405SHong Zhang 
163934ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);
1640a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(23));CHKERRQ(ierr);
164134ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1642b34f08ffSHong Zhang 
1643b34f08ffSHong Zhang       if (mumps->ninfo && mumps->ninfo <= 40){
1644b34f08ffSHong Zhang         PetscInt i;
1645b34f08ffSHong Zhang         for (i=0; i<mumps->ninfo; i++){
1646b34f08ffSHong Zhang           ierr = PetscViewerASCIIPrintf(viewer, "  INFO(%d): \n",mumps->info[i]);CHKERRQ(ierr);
1647b34f08ffSHong Zhang           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(mumps->info[i]));CHKERRQ(ierr);
1648b34f08ffSHong Zhang           ierr = PetscViewerFlush(viewer);
1649b34f08ffSHong Zhang         }
1650b34f08ffSHong Zhang       }
1651b34f08ffSHong Zhang 
1652b34f08ffSHong Zhang 
16531575c14dSBarry Smith       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
1654f6c57405SHong Zhang 
1655a5e57a09SHong Zhang       if (!mumps->myid) { /* information from the host */
1656a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",mumps->id.RINFOG(1));CHKERRQ(ierr);
1657a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",mumps->id.RINFOG(2));CHKERRQ(ierr);
1658a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",mumps->id.RINFOG(3));CHKERRQ(ierr);
1659a5e57a09SHong 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);
1660f6c57405SHong Zhang 
1661a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(3));CHKERRQ(ierr);
1662a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(4));CHKERRQ(ierr);
1663a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",mumps->id.INFOG(5));CHKERRQ(ierr);
1664a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",mumps->id.INFOG(6));CHKERRQ(ierr);
1665a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively use after analysis): %d \n",mumps->id.INFOG(7));CHKERRQ(ierr);
1666a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",mumps->id.INFOG(8));CHKERRQ(ierr);
1667a5e57a09SHong 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);
1668a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",mumps->id.INFOG(10));CHKERRQ(ierr);
1669a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",mumps->id.INFOG(11));CHKERRQ(ierr);
1670a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",mumps->id.INFOG(12));CHKERRQ(ierr);
1671a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",mumps->id.INFOG(13));CHKERRQ(ierr);
1672a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",mumps->id.INFOG(14));CHKERRQ(ierr);
1673a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",mumps->id.INFOG(15));CHKERRQ(ierr);
1674a5e57a09SHong 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);
1675a5e57a09SHong 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);
1676a5e57a09SHong 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);
1677a5e57a09SHong 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);
1678a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",mumps->id.INFOG(20));CHKERRQ(ierr);
1679a5e57a09SHong 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);
1680a5e57a09SHong 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);
1681a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",mumps->id.INFOG(23));CHKERRQ(ierr);
1682a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",mumps->id.INFOG(24));CHKERRQ(ierr);
1683a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",mumps->id.INFOG(25));CHKERRQ(ierr);
168440d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(28) (after factorization: number of null pivots encountered): %d\n",mumps->id.INFOG(28));CHKERRQ(ierr);
168540d435e3SHong 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);
168640d435e3SHong 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);
168740d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(32) (after analysis: type of analysis done): %d\n",mumps->id.INFOG(32));CHKERRQ(ierr);
168840d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(33) (value used for ICNTL(8)): %d\n",mumps->id.INFOG(33));CHKERRQ(ierr);
168940d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(34) (exponent of the determinant if determinant is requested): %d\n",mumps->id.INFOG(34));CHKERRQ(ierr);
1690f6c57405SHong Zhang       }
1691f6c57405SHong Zhang     }
1692cb828f0fSHong Zhang   }
1693f6c57405SHong Zhang   PetscFunctionReturn(0);
1694f6c57405SHong Zhang }
1695f6c57405SHong Zhang 
169635bd34faSBarry Smith #undef __FUNCT__
169735bd34faSBarry Smith #define __FUNCT__ "MatGetInfo_MUMPS"
169835bd34faSBarry Smith PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
169935bd34faSBarry Smith {
1700cb828f0fSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)A->spptr;
170135bd34faSBarry Smith 
170235bd34faSBarry Smith   PetscFunctionBegin;
170335bd34faSBarry Smith   info->block_size        = 1.0;
1704cb828f0fSHong Zhang   info->nz_allocated      = mumps->id.INFOG(20);
1705cb828f0fSHong Zhang   info->nz_used           = mumps->id.INFOG(20);
170635bd34faSBarry Smith   info->nz_unneeded       = 0.0;
170735bd34faSBarry Smith   info->assemblies        = 0.0;
170835bd34faSBarry Smith   info->mallocs           = 0.0;
170935bd34faSBarry Smith   info->memory            = 0.0;
171035bd34faSBarry Smith   info->fill_ratio_given  = 0;
171135bd34faSBarry Smith   info->fill_ratio_needed = 0;
171235bd34faSBarry Smith   info->factor_mallocs    = 0;
171335bd34faSBarry Smith   PetscFunctionReturn(0);
171435bd34faSBarry Smith }
171535bd34faSBarry Smith 
17165ccb76cbSHong Zhang /* -------------------------------------------------------------------------------------------*/
17175ccb76cbSHong Zhang #undef __FUNCT__
17188e7ba810SStefano Zampini #define __FUNCT__ "MatFactorSetSchurIS_MUMPS"
17198e7ba810SStefano Zampini PetscErrorCode MatFactorSetSchurIS_MUMPS(Mat F, IS is)
17206444a565SStefano Zampini {
17216444a565SStefano Zampini   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
17228e7ba810SStefano Zampini   const PetscInt *idxs;
17238e7ba810SStefano Zampini   PetscInt       size,i;
17246444a565SStefano Zampini   PetscErrorCode ierr;
17256444a565SStefano Zampini 
17266444a565SStefano Zampini   PetscFunctionBegin;
172759ac8732SStefano Zampini   if (mumps->size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MUMPS parallel Schur complements not yet supported from PETSc\n");
17288e7ba810SStefano Zampini   ierr = ISGetLocalSize(is,&size);CHKERRQ(ierr);
17296444a565SStefano Zampini   if (mumps->id.size_schur != size) {
17306444a565SStefano Zampini     ierr = PetscFree2(mumps->id.listvar_schur,mumps->id.schur);CHKERRQ(ierr);
17316444a565SStefano Zampini     mumps->id.size_schur = size;
17326444a565SStefano Zampini     mumps->id.schur_lld = size;
17336444a565SStefano Zampini     ierr = PetscMalloc2(size,&mumps->id.listvar_schur,size*size,&mumps->id.schur);CHKERRQ(ierr);
17346444a565SStefano Zampini   }
17358e7ba810SStefano Zampini   ierr = ISGetIndices(is,&idxs);CHKERRQ(ierr);
17366444a565SStefano Zampini   ierr = PetscMemcpy(mumps->id.listvar_schur,idxs,size*sizeof(PetscInt));CHKERRQ(ierr);
17378e7ba810SStefano Zampini   /* MUMPS expects Fortran style indices */
17388e7ba810SStefano Zampini   for (i=0;i<size;i++) mumps->id.listvar_schur[i]++;
17398e7ba810SStefano Zampini   ierr = ISRestoreIndices(is,&idxs);CHKERRQ(ierr);
174078a7176cSStefano Zampini   if (size) { /* turn on Schur switch if we the set of indices is not empty */
17416444a565SStefano Zampini     if (F->factortype == MAT_FACTOR_LU) {
174259ac8732SStefano Zampini       mumps->id.ICNTL(19) = 3; /* MUMPS returns full matrix */
17436444a565SStefano Zampini     } else {
174459ac8732SStefano Zampini       mumps->id.ICNTL(19) = 2; /* MUMPS returns lower triangular part */
17456444a565SStefano Zampini     }
174659ac8732SStefano Zampini     /* set a special value of ICNTL (not handled my MUMPS) to be used in the solve phase by PETSc */
1747b5fa320bSStefano Zampini     mumps->id.ICNTL(26) = -1;
17486444a565SStefano Zampini   }
17496444a565SStefano Zampini   PetscFunctionReturn(0);
17506444a565SStefano Zampini }
175159ac8732SStefano Zampini 
17526444a565SStefano Zampini /* -------------------------------------------------------------------------------------------*/
17536444a565SStefano Zampini #undef __FUNCT__
17545a05ddb0SStefano Zampini #define __FUNCT__ "MatFactorCreateSchurComplement_MUMPS"
17555a05ddb0SStefano Zampini PetscErrorCode MatFactorCreateSchurComplement_MUMPS(Mat F,Mat* S)
17566444a565SStefano Zampini {
17576444a565SStefano Zampini   Mat            St;
17586444a565SStefano Zampini   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
17596444a565SStefano Zampini   PetscScalar    *array;
17606444a565SStefano Zampini #if defined(PETSC_USE_COMPLEX)
17618ac429a0SStefano Zampini   PetscScalar    im = PetscSqrtScalar((PetscScalar)-1.0);
17626444a565SStefano Zampini #endif
17636444a565SStefano Zampini   PetscErrorCode ierr;
17646444a565SStefano Zampini 
17656444a565SStefano Zampini   PetscFunctionBegin;
17665a05ddb0SStefano 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");
17675a05ddb0SStefano Zampini   else if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before");
17686f3cc6f9SBarry Smith 
17696444a565SStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)F),&St);CHKERRQ(ierr);
17706444a565SStefano Zampini   ierr = MatSetSizes(St,PETSC_DECIDE,PETSC_DECIDE,mumps->id.size_schur,mumps->id.size_schur);CHKERRQ(ierr);
17716444a565SStefano Zampini   ierr = MatSetType(St,MATDENSE);CHKERRQ(ierr);
17726444a565SStefano Zampini   ierr = MatSetUp(St);CHKERRQ(ierr);
17736444a565SStefano Zampini   ierr = MatDenseGetArray(St,&array);CHKERRQ(ierr);
177459ac8732SStefano Zampini   if (!mumps->sym) { /* MUMPS always return a full matrix */
17756444a565SStefano Zampini     if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
17766444a565SStefano Zampini       PetscInt i,j,N=mumps->id.size_schur;
17776444a565SStefano Zampini       for (i=0;i<N;i++) {
17786444a565SStefano Zampini         for (j=0;j<N;j++) {
17796444a565SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
17806444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j];
17816444a565SStefano Zampini #else
17826444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
17836444a565SStefano Zampini #endif
17846444a565SStefano Zampini           array[j*N+i] = val;
17856444a565SStefano Zampini         }
17866444a565SStefano Zampini       }
17876444a565SStefano Zampini     } else { /* stored by columns */
17886444a565SStefano Zampini       ierr = PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
17896444a565SStefano Zampini     }
17906444a565SStefano Zampini   } else { /* either full or lower-triangular (not packed) */
17916444a565SStefano Zampini     if (mumps->id.ICNTL(19) == 2) { /* lower triangular stored by columns */
17926444a565SStefano Zampini       PetscInt i,j,N=mumps->id.size_schur;
17936444a565SStefano Zampini       for (i=0;i<N;i++) {
17946444a565SStefano Zampini         for (j=i;j<N;j++) {
17956444a565SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
17966444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j];
17976444a565SStefano Zampini #else
17986444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
17996444a565SStefano Zampini #endif
18006444a565SStefano Zampini           array[i*N+j] = val;
18016444a565SStefano Zampini           array[j*N+i] = val;
18026444a565SStefano Zampini         }
18036444a565SStefano Zampini       }
18046444a565SStefano Zampini     } else if (mumps->id.ICNTL(19) == 3) { /* full matrix */
18056444a565SStefano Zampini       ierr = PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
18066444a565SStefano Zampini     } else { /* ICNTL(19) == 1 lower triangular stored by rows */
18076444a565SStefano Zampini       PetscInt i,j,N=mumps->id.size_schur;
18086444a565SStefano Zampini       for (i=0;i<N;i++) {
18096444a565SStefano Zampini         for (j=0;j<i+1;j++) {
18106444a565SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
18116444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j];
18126444a565SStefano Zampini #else
18136444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
18146444a565SStefano Zampini #endif
18156444a565SStefano Zampini           array[i*N+j] = val;
18166444a565SStefano Zampini           array[j*N+i] = val;
18176444a565SStefano Zampini         }
18186444a565SStefano Zampini       }
18196444a565SStefano Zampini     }
18206444a565SStefano Zampini   }
18216444a565SStefano Zampini   ierr = MatDenseRestoreArray(St,&array);CHKERRQ(ierr);
18226444a565SStefano Zampini   *S = St;
18236444a565SStefano Zampini   PetscFunctionReturn(0);
18246444a565SStefano Zampini }
18256444a565SStefano Zampini 
182659ac8732SStefano Zampini /* -------------------------------------------------------------------------------------------*/
182759ac8732SStefano Zampini #undef __FUNCT__
18285a05ddb0SStefano Zampini #define __FUNCT__ "MatFactorGetSchurComplement_MUMPS"
18295a05ddb0SStefano Zampini PetscErrorCode MatFactorGetSchurComplement_MUMPS(Mat F,Mat* S)
183059ac8732SStefano Zampini {
183159ac8732SStefano Zampini   Mat            St;
183259ac8732SStefano Zampini   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
183359ac8732SStefano Zampini   PetscErrorCode ierr;
183459ac8732SStefano Zampini 
183559ac8732SStefano Zampini   PetscFunctionBegin;
18365a05ddb0SStefano 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");
18375a05ddb0SStefano Zampini   else if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before");
18386f3cc6f9SBarry Smith 
18397bc32933SStefano 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 */
184059ac8732SStefano Zampini   ierr = MatCreateSeqDense(PetscObjectComm((PetscObject)F),mumps->id.size_schur,mumps->id.size_schur,(PetscScalar*)mumps->id.schur,&St);CHKERRQ(ierr);
184159ac8732SStefano Zampini   *S = St;
184259ac8732SStefano Zampini   PetscFunctionReturn(0);
184359ac8732SStefano Zampini }
184459ac8732SStefano Zampini 
184559ac8732SStefano Zampini /* -------------------------------------------------------------------------------------------*/
184659ac8732SStefano Zampini #undef __FUNCT__
1847*6dba178dSStefano Zampini #define __FUNCT__ "MatFactorFactorizeSchurComplement_MUMPS"
1848*6dba178dSStefano Zampini PetscErrorCode MatFactorFactorizeSchurComplement_MUMPS(Mat F)
18493280fbe9SStefano Zampini {
18503280fbe9SStefano Zampini   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
18513280fbe9SStefano Zampini   PetscErrorCode ierr;
18523280fbe9SStefano Zampini 
18533280fbe9SStefano Zampini   PetscFunctionBegin;
18543280fbe9SStefano Zampini   if (!mumps->id.ICNTL(19)) { /* do nothing */
18553280fbe9SStefano Zampini     PetscFunctionReturn(0);
18563280fbe9SStefano Zampini   }
185734588265SStefano Zampini   if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatMumpsSetSchurIndices before");
18583280fbe9SStefano Zampini   ierr = MatMumpsFactorSchur_Private(mumps);CHKERRQ(ierr);
18593280fbe9SStefano Zampini   PetscFunctionReturn(0);
18603280fbe9SStefano Zampini }
18613280fbe9SStefano Zampini 
18623280fbe9SStefano Zampini #undef __FUNCT__
18635a05ddb0SStefano Zampini #define __FUNCT__ "MatFactorInvertSchurComplement_MUMPS"
18645a05ddb0SStefano Zampini PetscErrorCode MatFactorInvertSchurComplement_MUMPS(Mat F)
186559ac8732SStefano Zampini {
186659ac8732SStefano Zampini   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
186759ac8732SStefano Zampini   PetscErrorCode ierr;
186859ac8732SStefano Zampini 
186959ac8732SStefano Zampini   PetscFunctionBegin;
187059ac8732SStefano Zampini   if (!mumps->id.ICNTL(19)) { /* do nothing */
187159ac8732SStefano Zampini     PetscFunctionReturn(0);
187259ac8732SStefano Zampini   }
18735a05ddb0SStefano Zampini   if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before");
187459ac8732SStefano Zampini   ierr = MatMumpsInvertSchur_Private(mumps);CHKERRQ(ierr);
187559ac8732SStefano Zampini   PetscFunctionReturn(0);
187659ac8732SStefano Zampini }
187759ac8732SStefano Zampini 
18786444a565SStefano Zampini /* -------------------------------------------------------------------------------------------*/
18796444a565SStefano Zampini #undef __FUNCT__
18805a05ddb0SStefano Zampini #define __FUNCT__ "MatFactorSolveSchurComplement_MUMPS"
18815a05ddb0SStefano Zampini PetscErrorCode MatFactorSolveSchurComplement_MUMPS(Mat F, Vec rhs, Vec sol)
1882e807eca7SStefano Zampini {
1883e807eca7SStefano Zampini   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
1884e807eca7SStefano Zampini   MumpsScalar    *orhs;
1885e807eca7SStefano Zampini   PetscScalar    *osol,*nrhs,*nsol;
18869a3a5937SStefano Zampini   PetscInt       orhs_size,osol_size,olrhs_size;
1887e807eca7SStefano Zampini   PetscErrorCode ierr;
1888e807eca7SStefano Zampini 
1889e807eca7SStefano Zampini   PetscFunctionBegin;
18905a05ddb0SStefano 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");
18915a05ddb0SStefano Zampini   if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before");
18926f3cc6f9SBarry Smith 
1893e807eca7SStefano Zampini   /* swap pointers */
1894e807eca7SStefano Zampini   orhs = mumps->id.redrhs;
18959a3a5937SStefano Zampini   olrhs_size = mumps->id.lredrhs;
1896a12f35bfSStefano Zampini   orhs_size = mumps->sizeredrhs;
1897e807eca7SStefano Zampini   osol = mumps->schur_sol;
1898a12f35bfSStefano Zampini   osol_size = mumps->schur_sizesol;
1899e807eca7SStefano Zampini   ierr = VecGetArray(rhs,&nrhs);CHKERRQ(ierr);
1900e807eca7SStefano Zampini   ierr = VecGetArray(sol,&nsol);CHKERRQ(ierr);
1901e807eca7SStefano Zampini   mumps->id.redrhs = (MumpsScalar*)nrhs;
1902a12f35bfSStefano Zampini   ierr = VecGetLocalSize(rhs,&mumps->sizeredrhs);CHKERRQ(ierr);
19039a3a5937SStefano Zampini   mumps->id.lredrhs = mumps->sizeredrhs;
1904e807eca7SStefano Zampini   mumps->schur_sol = nsol;
1905a12f35bfSStefano Zampini   ierr = VecGetLocalSize(sol,&mumps->schur_sizesol);CHKERRQ(ierr);
1906a12f35bfSStefano Zampini 
1907e807eca7SStefano Zampini   /* solve Schur complement */
1908e807eca7SStefano Zampini   mumps->id.nrhs = 1;
1909e807eca7SStefano Zampini   ierr = MatMumpsSolveSchur_Private(mumps,PETSC_FALSE);CHKERRQ(ierr);
1910e807eca7SStefano Zampini   /* restore pointers */
1911e807eca7SStefano Zampini   ierr = VecRestoreArray(rhs,&nrhs);CHKERRQ(ierr);
1912e807eca7SStefano Zampini   ierr = VecRestoreArray(sol,&nsol);CHKERRQ(ierr);
1913e807eca7SStefano Zampini   mumps->id.redrhs = orhs;
19149a3a5937SStefano Zampini   mumps->id.lredrhs = olrhs_size;
1915a12f35bfSStefano Zampini   mumps->sizeredrhs = orhs_size;
1916e807eca7SStefano Zampini   mumps->schur_sol = osol;
1917a12f35bfSStefano Zampini   mumps->schur_sizesol = osol_size;
1918e807eca7SStefano Zampini   PetscFunctionReturn(0);
1919e807eca7SStefano Zampini }
1920e807eca7SStefano Zampini 
19217404bcfbSStefano Zampini /* -------------------------------------------------------------------------------------------*/
19227404bcfbSStefano Zampini #undef __FUNCT__
19235a05ddb0SStefano Zampini #define __FUNCT__ "MatFactorSolveSchurComplementTranspose_MUMPS"
19245a05ddb0SStefano Zampini PetscErrorCode MatFactorSolveSchurComplementTranspose_MUMPS(Mat F, Vec rhs, Vec sol)
19257404bcfbSStefano Zampini {
19267404bcfbSStefano Zampini   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
19277404bcfbSStefano Zampini   MumpsScalar    *orhs;
19287404bcfbSStefano Zampini   PetscScalar    *osol,*nrhs,*nsol;
1929a12f35bfSStefano Zampini   PetscInt       orhs_size,osol_size;
19307404bcfbSStefano Zampini   PetscErrorCode ierr;
19317404bcfbSStefano Zampini 
19327404bcfbSStefano Zampini   PetscFunctionBegin;
19335a05ddb0SStefano 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");
19345a05ddb0SStefano Zampini   else if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before");
19356f3cc6f9SBarry Smith 
19367404bcfbSStefano Zampini   /* swap pointers */
19377404bcfbSStefano Zampini   orhs = mumps->id.redrhs;
1938a12f35bfSStefano Zampini   orhs_size = mumps->sizeredrhs;
19397404bcfbSStefano Zampini   osol = mumps->schur_sol;
1940a12f35bfSStefano Zampini   osol_size = mumps->schur_sizesol;
19417404bcfbSStefano Zampini   ierr = VecGetArray(rhs,&nrhs);CHKERRQ(ierr);
19427404bcfbSStefano Zampini   ierr = VecGetArray(sol,&nsol);CHKERRQ(ierr);
19437404bcfbSStefano Zampini   mumps->id.redrhs = (MumpsScalar*)nrhs;
1944a12f35bfSStefano Zampini   ierr = VecGetLocalSize(rhs,&mumps->sizeredrhs);CHKERRQ(ierr);
19457404bcfbSStefano Zampini   mumps->schur_sol = nsol;
1946a12f35bfSStefano Zampini   ierr = VecGetLocalSize(sol,&mumps->schur_sizesol);CHKERRQ(ierr);
1947a12f35bfSStefano Zampini 
19487404bcfbSStefano Zampini   /* solve Schur complement */
19497404bcfbSStefano Zampini   mumps->id.nrhs = 1;
19507404bcfbSStefano Zampini   mumps->id.ICNTL(9) = 0;
19517404bcfbSStefano Zampini   ierr = MatMumpsSolveSchur_Private(mumps,PETSC_FALSE);CHKERRQ(ierr);
19527404bcfbSStefano Zampini   mumps->id.ICNTL(9) = 1;
19537404bcfbSStefano Zampini   /* restore pointers */
19547404bcfbSStefano Zampini   ierr = VecRestoreArray(rhs,&nrhs);CHKERRQ(ierr);
19557404bcfbSStefano Zampini   ierr = VecRestoreArray(sol,&nsol);CHKERRQ(ierr);
19567404bcfbSStefano Zampini   mumps->id.redrhs = orhs;
1957a12f35bfSStefano Zampini   mumps->sizeredrhs = orhs_size;
19587404bcfbSStefano Zampini   mumps->schur_sol = osol;
1959a12f35bfSStefano Zampini   mumps->schur_sizesol = osol_size;
19607404bcfbSStefano Zampini   PetscFunctionReturn(0);
19617404bcfbSStefano Zampini }
19627404bcfbSStefano Zampini 
19637404bcfbSStefano Zampini #undef __FUNCT__
1964a0b0af32SStefano Zampini #define __FUNCT__ "MatMumpsSchurComplementSetSym"
1965a0b0af32SStefano Zampini /*@
1966a0b0af32SStefano Zampini   MatMumpsSchurComplementSetSym - Set symmetric info for Schur complement
1967a0b0af32SStefano Zampini 
1968a0b0af32SStefano Zampini    Logically Collective on Mat
1969a0b0af32SStefano Zampini 
1970a0b0af32SStefano Zampini    Input Parameters:
1971a0b0af32SStefano Zampini +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1972a0b0af32SStefano Zampini -  sym - either 0 (non-symmetric), 1 (symmetric positive definite) or 2 (symmetric indefinite) following MUMPS convention
1973a0b0af32SStefano Zampini 
1974a0b0af32SStefano Zampini    Notes:
1975a0b0af32SStefano Zampini    The parameter is used to compute the correct factorization of the Schur complement matrices
1976a0b0af32SStefano Zampini    This could be useful in case the nature of the Schur complement is different from that of the matrix to be factored
1977a0b0af32SStefano Zampini    MUMPS Schur complement mode is currently implemented for sequential matrices.
1978a0b0af32SStefano Zampini 
1979a0b0af32SStefano Zampini    Level: advanced
1980a0b0af32SStefano Zampini 
1981a0b0af32SStefano Zampini    References: MUMPS Users' Guide
1982a0b0af32SStefano Zampini 
1983a0b0af32SStefano Zampini .seealso: MatGetFactor(), MatMumpsSetSchurIndices()
1984a0b0af32SStefano Zampini @*/
1985a0b0af32SStefano Zampini PetscErrorCode MatMumpsSchurComplementSetSym(Mat F, PetscInt sym)
1986a0b0af32SStefano Zampini {
1987a0b0af32SStefano Zampini   PetscErrorCode ierr;
1988a0b0af32SStefano Zampini 
1989a0b0af32SStefano Zampini   PetscFunctionBegin;
1990a0b0af32SStefano Zampini   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
1991a0b0af32SStefano Zampini   PetscValidLogicalCollectiveInt(F,sym,2);
1992a0b0af32SStefano Zampini   ierr = PetscUseMethod(F,"MatMumpsSchurComplementSetSym_C",(Mat,PetscInt),(F,sym));CHKERRQ(ierr);
1993a0b0af32SStefano Zampini   PetscFunctionReturn(0);
1994a0b0af32SStefano Zampini }
1995a0b0af32SStefano Zampini 
1996a0b0af32SStefano Zampini /* -------------------------------------------------------------------------------------------*/
1997a0b0af32SStefano Zampini #undef __FUNCT__
1998a0b0af32SStefano Zampini #define __FUNCT__ "MatMumpsSchurComplementSetSym_MUMPS"
1999a0b0af32SStefano Zampini PetscErrorCode MatMumpsSchurComplementSetSym_MUMPS(Mat F, PetscInt sym)
2000a0b0af32SStefano Zampini {
2001a0b0af32SStefano Zampini   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
2002a0b0af32SStefano Zampini 
2003a0b0af32SStefano Zampini   PetscFunctionBegin;
2004a0b0af32SStefano Zampini   if (mumps->schur_factored && mumps->sym != mumps->schur_sym) {
2005a0b0af32SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONG,"Schur complement data has been already factored");
2006a0b0af32SStefano Zampini   }
2007a0b0af32SStefano Zampini   mumps->schur_sym = sym;
2008a0b0af32SStefano Zampini   PetscFunctionReturn(0);
2009a0b0af32SStefano Zampini }
2010a0b0af32SStefano Zampini 
2011e807eca7SStefano Zampini /* -------------------------------------------------------------------------------------------*/
2012e807eca7SStefano Zampini #undef __FUNCT__
20135ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl_MUMPS"
20145ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival)
20155ccb76cbSHong Zhang {
2016a5e57a09SHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
20175ccb76cbSHong Zhang 
20185ccb76cbSHong Zhang   PetscFunctionBegin;
2019a5e57a09SHong Zhang   mumps->id.ICNTL(icntl) = ival;
20205ccb76cbSHong Zhang   PetscFunctionReturn(0);
20215ccb76cbSHong Zhang }
20225ccb76cbSHong Zhang 
20235ccb76cbSHong Zhang #undef __FUNCT__
2024bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetIcntl_MUMPS"
2025bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt *ival)
2026bc6112feSHong Zhang {
2027bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
2028bc6112feSHong Zhang 
2029bc6112feSHong Zhang   PetscFunctionBegin;
2030bc6112feSHong Zhang   *ival = mumps->id.ICNTL(icntl);
2031bc6112feSHong Zhang   PetscFunctionReturn(0);
2032bc6112feSHong Zhang }
2033bc6112feSHong Zhang 
2034bc6112feSHong Zhang #undef __FUNCT__
20355ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl"
20365ccb76cbSHong Zhang /*@
20375ccb76cbSHong Zhang   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()
20385ccb76cbSHong Zhang 
20395ccb76cbSHong Zhang    Logically Collective on Mat
20405ccb76cbSHong Zhang 
20415ccb76cbSHong Zhang    Input Parameters:
20425ccb76cbSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
20435ccb76cbSHong Zhang .  icntl - index of MUMPS parameter array ICNTL()
20445ccb76cbSHong Zhang -  ival - value of MUMPS ICNTL(icntl)
20455ccb76cbSHong Zhang 
20465ccb76cbSHong Zhang   Options Database:
20475ccb76cbSHong Zhang .   -mat_mumps_icntl_<icntl> <ival>
20485ccb76cbSHong Zhang 
20495ccb76cbSHong Zhang    Level: beginner
20505ccb76cbSHong Zhang 
20515ccb76cbSHong Zhang    References: MUMPS Users' Guide
20525ccb76cbSHong Zhang 
20535ccb76cbSHong Zhang .seealso: MatGetFactor()
20545ccb76cbSHong Zhang @*/
20555ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival)
20565ccb76cbSHong Zhang {
20575ccb76cbSHong Zhang   PetscErrorCode ierr;
20585ccb76cbSHong Zhang 
20595ccb76cbSHong Zhang   PetscFunctionBegin;
20605ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
20615ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,ival,3);
20625ccb76cbSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr);
20635ccb76cbSHong Zhang   PetscFunctionReturn(0);
20645ccb76cbSHong Zhang }
20655ccb76cbSHong Zhang 
2066bc6112feSHong Zhang #undef __FUNCT__
2067bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetIcntl"
2068a21f80fcSHong Zhang /*@
2069a21f80fcSHong Zhang   MatMumpsGetIcntl - Get MUMPS parameter ICNTL()
2070a21f80fcSHong Zhang 
2071a21f80fcSHong Zhang    Logically Collective on Mat
2072a21f80fcSHong Zhang 
2073a21f80fcSHong Zhang    Input Parameters:
2074a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2075a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array ICNTL()
2076a21f80fcSHong Zhang 
2077a21f80fcSHong Zhang   Output Parameter:
2078a21f80fcSHong Zhang .  ival - value of MUMPS ICNTL(icntl)
2079a21f80fcSHong Zhang 
2080a21f80fcSHong Zhang    Level: beginner
2081a21f80fcSHong Zhang 
2082a21f80fcSHong Zhang    References: MUMPS Users' Guide
2083a21f80fcSHong Zhang 
2084a21f80fcSHong Zhang .seealso: MatGetFactor()
2085a21f80fcSHong Zhang @*/
2086bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl(Mat F,PetscInt icntl,PetscInt *ival)
2087bc6112feSHong Zhang {
2088bc6112feSHong Zhang   PetscErrorCode ierr;
2089bc6112feSHong Zhang 
2090bc6112feSHong Zhang   PetscFunctionBegin;
2091bc6112feSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
2092bc6112feSHong Zhang   PetscValidIntPointer(ival,3);
2093bc6112feSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetIcntl_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
2094bc6112feSHong Zhang   PetscFunctionReturn(0);
2095bc6112feSHong Zhang }
2096bc6112feSHong Zhang 
20978928b65cSHong Zhang /* -------------------------------------------------------------------------------------------*/
20988928b65cSHong Zhang #undef __FUNCT__
20998928b65cSHong Zhang #define __FUNCT__ "MatMumpsSetCntl_MUMPS"
21008928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal val)
21018928b65cSHong Zhang {
21028928b65cSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
21038928b65cSHong Zhang 
21048928b65cSHong Zhang   PetscFunctionBegin;
21058928b65cSHong Zhang   mumps->id.CNTL(icntl) = val;
21068928b65cSHong Zhang   PetscFunctionReturn(0);
21078928b65cSHong Zhang }
21088928b65cSHong Zhang 
21098928b65cSHong Zhang #undef __FUNCT__
2110bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetCntl_MUMPS"
2111bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal *val)
2112bc6112feSHong Zhang {
2113bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
2114bc6112feSHong Zhang 
2115bc6112feSHong Zhang   PetscFunctionBegin;
2116bc6112feSHong Zhang   *val = mumps->id.CNTL(icntl);
2117bc6112feSHong Zhang   PetscFunctionReturn(0);
2118bc6112feSHong Zhang }
2119bc6112feSHong Zhang 
2120bc6112feSHong Zhang #undef __FUNCT__
21218928b65cSHong Zhang #define __FUNCT__ "MatMumpsSetCntl"
21228928b65cSHong Zhang /*@
21238928b65cSHong Zhang   MatMumpsSetCntl - Set MUMPS parameter CNTL()
21248928b65cSHong Zhang 
21258928b65cSHong Zhang    Logically Collective on Mat
21268928b65cSHong Zhang 
21278928b65cSHong Zhang    Input Parameters:
21288928b65cSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
21298928b65cSHong Zhang .  icntl - index of MUMPS parameter array CNTL()
21308928b65cSHong Zhang -  val - value of MUMPS CNTL(icntl)
21318928b65cSHong Zhang 
21328928b65cSHong Zhang   Options Database:
21338928b65cSHong Zhang .   -mat_mumps_cntl_<icntl> <val>
21348928b65cSHong Zhang 
21358928b65cSHong Zhang    Level: beginner
21368928b65cSHong Zhang 
21378928b65cSHong Zhang    References: MUMPS Users' Guide
21388928b65cSHong Zhang 
21398928b65cSHong Zhang .seealso: MatGetFactor()
21408928b65cSHong Zhang @*/
21418928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl(Mat F,PetscInt icntl,PetscReal val)
21428928b65cSHong Zhang {
21438928b65cSHong Zhang   PetscErrorCode ierr;
21448928b65cSHong Zhang 
21458928b65cSHong Zhang   PetscFunctionBegin;
21468928b65cSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
2147bc6112feSHong Zhang   PetscValidLogicalCollectiveReal(F,val,3);
21488928b65cSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsSetCntl_C",(Mat,PetscInt,PetscReal),(F,icntl,val));CHKERRQ(ierr);
21498928b65cSHong Zhang   PetscFunctionReturn(0);
21508928b65cSHong Zhang }
21518928b65cSHong Zhang 
2152bc6112feSHong Zhang #undef __FUNCT__
2153bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetCntl"
2154a21f80fcSHong Zhang /*@
2155a21f80fcSHong Zhang   MatMumpsGetCntl - Get MUMPS parameter CNTL()
2156a21f80fcSHong Zhang 
2157a21f80fcSHong Zhang    Logically Collective on Mat
2158a21f80fcSHong Zhang 
2159a21f80fcSHong Zhang    Input Parameters:
2160a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2161a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array CNTL()
2162a21f80fcSHong Zhang 
2163a21f80fcSHong Zhang   Output Parameter:
2164a21f80fcSHong Zhang .  val - value of MUMPS CNTL(icntl)
2165a21f80fcSHong Zhang 
2166a21f80fcSHong Zhang    Level: beginner
2167a21f80fcSHong Zhang 
2168a21f80fcSHong Zhang    References: MUMPS Users' Guide
2169a21f80fcSHong Zhang 
2170a21f80fcSHong Zhang .seealso: MatGetFactor()
2171a21f80fcSHong Zhang @*/
2172bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl(Mat F,PetscInt icntl,PetscReal *val)
2173bc6112feSHong Zhang {
2174bc6112feSHong Zhang   PetscErrorCode ierr;
2175bc6112feSHong Zhang 
2176bc6112feSHong Zhang   PetscFunctionBegin;
2177bc6112feSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
2178bc6112feSHong Zhang   PetscValidRealPointer(val,3);
2179bc6112feSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetCntl_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
2180bc6112feSHong Zhang   PetscFunctionReturn(0);
2181bc6112feSHong Zhang }
2182bc6112feSHong Zhang 
2183bc6112feSHong Zhang #undef __FUNCT__
2184ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfo_MUMPS"
2185ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F,PetscInt icntl,PetscInt *info)
2186bc6112feSHong Zhang {
2187bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
2188bc6112feSHong Zhang 
2189bc6112feSHong Zhang   PetscFunctionBegin;
2190bc6112feSHong Zhang   *info = mumps->id.INFO(icntl);
2191bc6112feSHong Zhang   PetscFunctionReturn(0);
2192bc6112feSHong Zhang }
2193bc6112feSHong Zhang 
2194bc6112feSHong Zhang #undef __FUNCT__
2195ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfog_MUMPS"
2196ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F,PetscInt icntl,PetscInt *infog)
2197bc6112feSHong Zhang {
2198bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
2199bc6112feSHong Zhang 
2200bc6112feSHong Zhang   PetscFunctionBegin;
2201bc6112feSHong Zhang   *infog = mumps->id.INFOG(icntl);
2202bc6112feSHong Zhang   PetscFunctionReturn(0);
2203bc6112feSHong Zhang }
2204bc6112feSHong Zhang 
2205bc6112feSHong Zhang #undef __FUNCT__
2206ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfo_MUMPS"
2207ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfo)
2208bc6112feSHong Zhang {
2209bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
2210bc6112feSHong Zhang 
2211bc6112feSHong Zhang   PetscFunctionBegin;
2212bc6112feSHong Zhang   *rinfo = mumps->id.RINFO(icntl);
2213bc6112feSHong Zhang   PetscFunctionReturn(0);
2214bc6112feSHong Zhang }
2215bc6112feSHong Zhang 
2216bc6112feSHong Zhang #undef __FUNCT__
2217ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfog_MUMPS"
2218ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfog)
2219bc6112feSHong Zhang {
2220bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
2221bc6112feSHong Zhang 
2222bc6112feSHong Zhang   PetscFunctionBegin;
2223bc6112feSHong Zhang   *rinfog = mumps->id.RINFOG(icntl);
2224bc6112feSHong Zhang   PetscFunctionReturn(0);
2225bc6112feSHong Zhang }
2226bc6112feSHong Zhang 
2227bc6112feSHong Zhang #undef __FUNCT__
2228ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfo"
2229a21f80fcSHong Zhang /*@
2230a21f80fcSHong Zhang   MatMumpsGetInfo - Get MUMPS parameter INFO()
2231a21f80fcSHong Zhang 
2232a21f80fcSHong Zhang    Logically Collective on Mat
2233a21f80fcSHong Zhang 
2234a21f80fcSHong Zhang    Input Parameters:
2235a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2236a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array INFO()
2237a21f80fcSHong Zhang 
2238a21f80fcSHong Zhang   Output Parameter:
2239a21f80fcSHong Zhang .  ival - value of MUMPS INFO(icntl)
2240a21f80fcSHong Zhang 
2241a21f80fcSHong Zhang    Level: beginner
2242a21f80fcSHong Zhang 
2243a21f80fcSHong Zhang    References: MUMPS Users' Guide
2244a21f80fcSHong Zhang 
2245a21f80fcSHong Zhang .seealso: MatGetFactor()
2246a21f80fcSHong Zhang @*/
2247ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo(Mat F,PetscInt icntl,PetscInt *ival)
2248bc6112feSHong Zhang {
2249bc6112feSHong Zhang   PetscErrorCode ierr;
2250bc6112feSHong Zhang 
2251bc6112feSHong Zhang   PetscFunctionBegin;
2252ca810319SHong Zhang   PetscValidIntPointer(ival,3);
2253ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetInfo_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
2254bc6112feSHong Zhang   PetscFunctionReturn(0);
2255bc6112feSHong Zhang }
2256bc6112feSHong Zhang 
2257bc6112feSHong Zhang #undef __FUNCT__
2258ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfog"
2259a21f80fcSHong Zhang /*@
2260a21f80fcSHong Zhang   MatMumpsGetInfog - Get MUMPS parameter INFOG()
2261a21f80fcSHong Zhang 
2262a21f80fcSHong Zhang    Logically Collective on Mat
2263a21f80fcSHong Zhang 
2264a21f80fcSHong Zhang    Input Parameters:
2265a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2266a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array INFOG()
2267a21f80fcSHong Zhang 
2268a21f80fcSHong Zhang   Output Parameter:
2269a21f80fcSHong Zhang .  ival - value of MUMPS INFOG(icntl)
2270a21f80fcSHong Zhang 
2271a21f80fcSHong Zhang    Level: beginner
2272a21f80fcSHong Zhang 
2273a21f80fcSHong Zhang    References: MUMPS Users' Guide
2274a21f80fcSHong Zhang 
2275a21f80fcSHong Zhang .seealso: MatGetFactor()
2276a21f80fcSHong Zhang @*/
2277ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog(Mat F,PetscInt icntl,PetscInt *ival)
2278bc6112feSHong Zhang {
2279bc6112feSHong Zhang   PetscErrorCode ierr;
2280bc6112feSHong Zhang 
2281bc6112feSHong Zhang   PetscFunctionBegin;
2282ca810319SHong Zhang   PetscValidIntPointer(ival,3);
2283ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetInfog_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
2284bc6112feSHong Zhang   PetscFunctionReturn(0);
2285bc6112feSHong Zhang }
2286bc6112feSHong Zhang 
2287bc6112feSHong Zhang #undef __FUNCT__
2288ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfo"
2289a21f80fcSHong Zhang /*@
2290a21f80fcSHong Zhang   MatMumpsGetRinfo - Get MUMPS parameter RINFO()
2291a21f80fcSHong Zhang 
2292a21f80fcSHong Zhang    Logically Collective on Mat
2293a21f80fcSHong Zhang 
2294a21f80fcSHong Zhang    Input Parameters:
2295a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2296a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array RINFO()
2297a21f80fcSHong Zhang 
2298a21f80fcSHong Zhang   Output Parameter:
2299a21f80fcSHong Zhang .  val - value of MUMPS RINFO(icntl)
2300a21f80fcSHong Zhang 
2301a21f80fcSHong Zhang    Level: beginner
2302a21f80fcSHong Zhang 
2303a21f80fcSHong Zhang    References: MUMPS Users' Guide
2304a21f80fcSHong Zhang 
2305a21f80fcSHong Zhang .seealso: MatGetFactor()
2306a21f80fcSHong Zhang @*/
2307ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo(Mat F,PetscInt icntl,PetscReal *val)
2308bc6112feSHong Zhang {
2309bc6112feSHong Zhang   PetscErrorCode ierr;
2310bc6112feSHong Zhang 
2311bc6112feSHong Zhang   PetscFunctionBegin;
2312bc6112feSHong Zhang   PetscValidRealPointer(val,3);
2313ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetRinfo_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
2314bc6112feSHong Zhang   PetscFunctionReturn(0);
2315bc6112feSHong Zhang }
2316bc6112feSHong Zhang 
2317bc6112feSHong Zhang #undef __FUNCT__
2318ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfog"
2319a21f80fcSHong Zhang /*@
2320a21f80fcSHong Zhang   MatMumpsGetRinfog - Get MUMPS parameter RINFOG()
2321a21f80fcSHong Zhang 
2322a21f80fcSHong Zhang    Logically Collective on Mat
2323a21f80fcSHong Zhang 
2324a21f80fcSHong Zhang    Input Parameters:
2325a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2326a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array RINFOG()
2327a21f80fcSHong Zhang 
2328a21f80fcSHong Zhang   Output Parameter:
2329a21f80fcSHong Zhang .  val - value of MUMPS RINFOG(icntl)
2330a21f80fcSHong Zhang 
2331a21f80fcSHong Zhang    Level: beginner
2332a21f80fcSHong Zhang 
2333a21f80fcSHong Zhang    References: MUMPS Users' Guide
2334a21f80fcSHong Zhang 
2335a21f80fcSHong Zhang .seealso: MatGetFactor()
2336a21f80fcSHong Zhang @*/
2337ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog(Mat F,PetscInt icntl,PetscReal *val)
2338bc6112feSHong Zhang {
2339bc6112feSHong Zhang   PetscErrorCode ierr;
2340bc6112feSHong Zhang 
2341bc6112feSHong Zhang   PetscFunctionBegin;
2342bc6112feSHong Zhang   PetscValidRealPointer(val,3);
2343ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetRinfog_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
2344bc6112feSHong Zhang   PetscFunctionReturn(0);
2345bc6112feSHong Zhang }
2346bc6112feSHong Zhang 
234724b6179bSKris Buschelman /*MC
23482692d6eeSBarry Smith   MATSOLVERMUMPS -  A matrix type providing direct solvers (LU and Cholesky) for
234924b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
235024b6179bSKris Buschelman 
235141c8de11SBarry Smith   Works with MATAIJ and MATSBAIJ matrices
235224b6179bSKris Buschelman 
2353c2b89b5dSBarry Smith   Use ./configure --download-mumps --download-scalapack --download-parmetis --download-metis --download-ptscotch  to have PETSc installed with MUMPS
2354c2b89b5dSBarry Smith 
2355c2b89b5dSBarry Smith   Use -pc_type cholesky or lu -pc_factor_mat_solver_package mumps to us this direct solver
2356c2b89b5dSBarry Smith 
235724b6179bSKris Buschelman   Options Database Keys:
23584e34a73bSHong Zhang +  -mat_mumps_icntl_1 <6>: ICNTL(1): output stream for error messages (None)
23594e34a73bSHong Zhang .  -mat_mumps_icntl_2 <0>: ICNTL(2): output stream for diagnostic printing, statistics, and warning (None)
23604e34a73bSHong Zhang .  -mat_mumps_icntl_3 <0>: ICNTL(3): output stream for global information, collected on the host (None)
23614e34a73bSHong Zhang .  -mat_mumps_icntl_4 <0>: ICNTL(4): level of printing (0 to 4) (None)
23624e34a73bSHong Zhang .  -mat_mumps_icntl_6 <7>: ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7) (None)
23634e34a73bSHong Zhang .  -mat_mumps_icntl_7 <7>: ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis (None)
23644e34a73bSHong Zhang .  -mat_mumps_icntl_8 <77>: ICNTL(8): scaling strategy (-2 to 8 or 77) (None)
23654e34a73bSHong Zhang .  -mat_mumps_icntl_10 <0>: ICNTL(10): max num of refinements (None)
23664e34a73bSHong Zhang .  -mat_mumps_icntl_11 <0>: ICNTL(11): statistics related to an error analysis (via -ksp_view) (None)
23674e34a73bSHong Zhang .  -mat_mumps_icntl_12 <1>: ICNTL(12): an ordering strategy for symmetric matrices (0 to 3) (None)
23684e34a73bSHong Zhang .  -mat_mumps_icntl_13 <0>: ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting (None)
23694e34a73bSHong Zhang .  -mat_mumps_icntl_14 <20>: ICNTL(14): percentage increase in the estimated working space (None)
23704e34a73bSHong Zhang .  -mat_mumps_icntl_19 <0>: ICNTL(19): computes the Schur complement (None)
23714e34a73bSHong Zhang .  -mat_mumps_icntl_22 <0>: ICNTL(22): in-core/out-of-core factorization and solve (0 or 1) (None)
23724e34a73bSHong Zhang .  -mat_mumps_icntl_23 <0>: ICNTL(23): max size of the working memory (MB) that can allocate per processor (None)
23734e34a73bSHong Zhang .  -mat_mumps_icntl_24 <0>: ICNTL(24): detection of null pivot rows (0 or 1) (None)
23744e34a73bSHong Zhang .  -mat_mumps_icntl_25 <0>: ICNTL(25): compute a solution of a deficient matrix and a null space basis (None)
23754e34a73bSHong Zhang .  -mat_mumps_icntl_26 <0>: ICNTL(26): drives the solution phase if a Schur complement matrix (None)
23764e34a73bSHong Zhang .  -mat_mumps_icntl_28 <1>: ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering (None)
23774e34a73bSHong Zhang .  -mat_mumps_icntl_29 <0>: ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis (None)
23784e34a73bSHong Zhang .  -mat_mumps_icntl_30 <0>: ICNTL(30): compute user-specified set of entries in inv(A) (None)
23794e34a73bSHong Zhang .  -mat_mumps_icntl_31 <0>: ICNTL(31): indicates which factors may be discarded during factorization (None)
23804e34a73bSHong Zhang .  -mat_mumps_icntl_33 <0>: ICNTL(33): compute determinant (None)
23814e34a73bSHong Zhang .  -mat_mumps_cntl_1 <0.01>: CNTL(1): relative pivoting threshold (None)
23824e34a73bSHong Zhang .  -mat_mumps_cntl_2 <1.49012e-08>: CNTL(2): stopping criterion of refinement (None)
23834e34a73bSHong Zhang .  -mat_mumps_cntl_3 <0>: CNTL(3): absolute pivoting threshold (None)
23844e34a73bSHong Zhang .  -mat_mumps_cntl_4 <-1>: CNTL(4): value for static pivoting (None)
23854e34a73bSHong Zhang -  -mat_mumps_cntl_5 <0>: CNTL(5): fixation for null pivots (None)
238624b6179bSKris Buschelman 
238724b6179bSKris Buschelman   Level: beginner
238824b6179bSKris Buschelman 
238941c8de11SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage
239041c8de11SBarry Smith 
239124b6179bSKris Buschelman M*/
239224b6179bSKris Buschelman 
239335bd34faSBarry Smith #undef __FUNCT__
239435bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_mumps"
2395f7a08781SBarry Smith static PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type)
239635bd34faSBarry Smith {
239735bd34faSBarry Smith   PetscFunctionBegin;
23982692d6eeSBarry Smith   *type = MATSOLVERMUMPS;
239935bd34faSBarry Smith   PetscFunctionReturn(0);
240035bd34faSBarry Smith }
240135bd34faSBarry Smith 
2402bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI AIJ matrices */
24032877fffaSHong Zhang #undef __FUNCT__
2404bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_aij_mumps"
24058cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
24062877fffaSHong Zhang {
24072877fffaSHong Zhang   Mat            B;
24082877fffaSHong Zhang   PetscErrorCode ierr;
24092877fffaSHong Zhang   Mat_MUMPS      *mumps;
2410ace3abfcSBarry Smith   PetscBool      isSeqAIJ;
24112877fffaSHong Zhang 
24122877fffaSHong Zhang   PetscFunctionBegin;
24132877fffaSHong Zhang   /* Create the factorization matrix */
2414251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr);
2415ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
24162877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
24172877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
2418bccb9932SShri Abhyankar   if (isSeqAIJ) {
24190298fd71SBarry Smith     ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr);
2420bccb9932SShri Abhyankar   } else {
24210298fd71SBarry Smith     ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr);
2422bccb9932SShri Abhyankar   }
24232877fffaSHong Zhang 
2424b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
24252205254eSKarl Rupp 
24262877fffaSHong Zhang   B->ops->view        = MatView_MUMPS;
242735bd34faSBarry Smith   B->ops->getinfo     = MatGetInfo_MUMPS;
242820be8e61SHong Zhang   B->ops->getdiagonal = MatGetDiagonal_MUMPS;
24292205254eSKarl Rupp 
2430bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
24315a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr);
24325a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorInvertSchurComplement_C",MatFactorInvertSchurComplement_MUMPS);CHKERRQ(ierr);
24335a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr);
24345a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSchurComplement_C",MatFactorGetSchurComplement_MUMPS);CHKERRQ(ierr);
24355a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplement_C",MatFactorSolveSchurComplement_MUMPS);CHKERRQ(ierr);
24365a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplementTranspose_C",MatFactorSolveSchurComplementTranspose_MUMPS);CHKERRQ(ierr);
2437*6dba178dSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorFactorizeSchurComplement_C",MatFactorFactorizeSchurComplement_MUMPS);CHKERRQ(ierr);
2438*6dba178dSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSchurComplementSetSym_C",MatMumpsSchurComplementSetSym_MUMPS);CHKERRQ(ierr);
2439bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
2440bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
2441bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
2442bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
2443ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
2444ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
2445ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
2446ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
24476444a565SStefano Zampini 
2448450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
2449450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
2450d5f3da31SBarry Smith     B->factortype            = MAT_FACTOR_LU;
2451bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
2452bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
2453746480a1SHong Zhang     mumps->sym = 0;
2454dcd589f8SShri Abhyankar   } else {
245567877ebaSShri Abhyankar     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2456450b117fSShri Abhyankar     B->factortype                  = MAT_FACTOR_CHOLESKY;
2457bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
2458bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
245959ac8732SStefano Zampini #if defined(PETSC_USE_COMPLEX)
246059ac8732SStefano Zampini     mumps->sym = 2;
246159ac8732SStefano Zampini #else
24626fdc2a6dSBarry Smith     if (A->spd_set && A->spd) mumps->sym = 1;
24636fdc2a6dSBarry Smith     else                      mumps->sym = 2;
246459ac8732SStefano Zampini #endif
2465450b117fSShri Abhyankar   }
24662877fffaSHong Zhang 
24672877fffaSHong Zhang   mumps->isAIJ    = PETSC_TRUE;
2468bf0cc555SLisandro Dalcin   mumps->Destroy  = B->ops->destroy;
24692877fffaSHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
24702877fffaSHong Zhang   B->spptr        = (void*)mumps;
24712205254eSKarl Rupp 
2472f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
2473746480a1SHong Zhang 
24742877fffaSHong Zhang   *F = B;
24752877fffaSHong Zhang   PetscFunctionReturn(0);
24762877fffaSHong Zhang }
24772877fffaSHong Zhang 
2478bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI SBAIJ matrices */
24792877fffaSHong Zhang #undef __FUNCT__
2480bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_sbaij_mumps"
24818cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
24822877fffaSHong Zhang {
24832877fffaSHong Zhang   Mat            B;
24842877fffaSHong Zhang   PetscErrorCode ierr;
24852877fffaSHong Zhang   Mat_MUMPS      *mumps;
2486ace3abfcSBarry Smith   PetscBool      isSeqSBAIJ;
24872877fffaSHong Zhang 
24882877fffaSHong Zhang   PetscFunctionBegin;
2489ce94432eSBarry Smith   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
2490ce94432eSBarry 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");
2491251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
24922877fffaSHong Zhang   /* Create the factorization matrix */
2493ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
24942877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
24952877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
2496b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
2497bccb9932SShri Abhyankar   if (isSeqSBAIJ) {
24980298fd71SBarry Smith     ierr = MatSeqSBAIJSetPreallocation(B,1,0,NULL);CHKERRQ(ierr);
24992205254eSKarl Rupp 
250016ebf90aSShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
2501dcd589f8SShri Abhyankar   } else {
25020298fd71SBarry Smith     ierr = MatMPISBAIJSetPreallocation(B,1,0,NULL,0,NULL);CHKERRQ(ierr);
25032205254eSKarl Rupp 
2504bccb9932SShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
2505bccb9932SShri Abhyankar   }
2506bccb9932SShri Abhyankar 
250767877ebaSShri Abhyankar   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2508bccb9932SShri Abhyankar   B->ops->view                   = MatView_MUMPS;
250920be8e61SHong Zhang   B->ops->getdiagonal            = MatGetDiagonal_MUMPS;
25102205254eSKarl Rupp 
2511bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
25125a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr);
25135a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorInvertSchurComplement_C",MatFactorInvertSchurComplement_MUMPS);CHKERRQ(ierr);
25145a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr);
25155a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSchurComplement_C",MatFactorGetSchurComplement_MUMPS);CHKERRQ(ierr);
25165a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplement_C",MatFactorSolveSchurComplement_MUMPS);CHKERRQ(ierr);
25175a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplementTranspose_C",MatFactorSolveSchurComplementTranspose_MUMPS);CHKERRQ(ierr);
2518*6dba178dSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorFactorizeSchurComplement_C",MatFactorFactorizeSchurComplement_MUMPS);CHKERRQ(ierr);
2519*6dba178dSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSchurComplementSetSym_C",MatMumpsSchurComplementSetSym_MUMPS);CHKERRQ(ierr);
2520b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
2521b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
2522b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
2523b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
2524ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
2525ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
2526ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
2527ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
25282205254eSKarl Rupp 
2529f4762488SHong Zhang   B->factortype = MAT_FACTOR_CHOLESKY;
253059ac8732SStefano Zampini #if defined(PETSC_USE_COMPLEX)
253159ac8732SStefano Zampini   mumps->sym = 2;
253259ac8732SStefano Zampini #else
25336fdc2a6dSBarry Smith   if (A->spd_set && A->spd) mumps->sym = 1;
25346fdc2a6dSBarry Smith   else                      mumps->sym = 2;
253559ac8732SStefano Zampini #endif
2536a214ac2aSShri Abhyankar 
2537bccb9932SShri Abhyankar   mumps->isAIJ    = PETSC_FALSE;
2538bf0cc555SLisandro Dalcin   mumps->Destroy  = B->ops->destroy;
2539f3c0ef26SHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
25402877fffaSHong Zhang   B->spptr        = (void*)mumps;
25412205254eSKarl Rupp 
2542f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
2543746480a1SHong Zhang 
25442877fffaSHong Zhang   *F = B;
25452877fffaSHong Zhang   PetscFunctionReturn(0);
25462877fffaSHong Zhang }
254797969023SHong Zhang 
2548450b117fSShri Abhyankar #undef __FUNCT__
2549bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_baij_mumps"
25508cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
255167877ebaSShri Abhyankar {
255267877ebaSShri Abhyankar   Mat            B;
255367877ebaSShri Abhyankar   PetscErrorCode ierr;
255467877ebaSShri Abhyankar   Mat_MUMPS      *mumps;
2555ace3abfcSBarry Smith   PetscBool      isSeqBAIJ;
255667877ebaSShri Abhyankar 
255767877ebaSShri Abhyankar   PetscFunctionBegin;
255867877ebaSShri Abhyankar   /* Create the factorization matrix */
2559251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr);
2560ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
256167877ebaSShri Abhyankar   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
256267877ebaSShri Abhyankar   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
2563bccb9932SShri Abhyankar   if (isSeqBAIJ) {
25640298fd71SBarry Smith     ierr = MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,NULL);CHKERRQ(ierr);
2565bccb9932SShri Abhyankar   } else {
25660298fd71SBarry Smith     ierr = MatMPIBAIJSetPreallocation(B,A->rmap->bs,0,NULL,0,NULL);CHKERRQ(ierr);
2567bccb9932SShri Abhyankar   }
2568450b117fSShri Abhyankar 
2569b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
2570450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
2571450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
2572450b117fSShri Abhyankar     B->factortype            = MAT_FACTOR_LU;
2573bccb9932SShri Abhyankar     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
2574bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
2575746480a1SHong Zhang     mumps->sym = 0;
2576f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n");
2577bccb9932SShri Abhyankar 
2578450b117fSShri Abhyankar   B->ops->view        = MatView_MUMPS;
257920be8e61SHong Zhang   B->ops->getdiagonal = MatGetDiagonal_MUMPS;
25802205254eSKarl Rupp 
2581bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
25825a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr);
25835a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorInvertSchurComplement_C",MatFactorInvertSchurComplement_MUMPS);CHKERRQ(ierr);
25845a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr);
25855a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSchurComplement_C",MatFactorGetSchurComplement_MUMPS);CHKERRQ(ierr);
25865a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplement_C",MatFactorSolveSchurComplement_MUMPS);CHKERRQ(ierr);
25875a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplementTranspose_C",MatFactorSolveSchurComplementTranspose_MUMPS);CHKERRQ(ierr);
2588*6dba178dSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorFactorizeSchurComplement_C",MatFactorFactorizeSchurComplement_MUMPS);CHKERRQ(ierr);
2589*6dba178dSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSchurComplementSetSym_C",MatMumpsSchurComplementSetSym_MUMPS);CHKERRQ(ierr);
2590bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
2591bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
2592bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
2593bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
2594ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
2595ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
2596ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
2597ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
2598450b117fSShri Abhyankar 
2599450b117fSShri Abhyankar   mumps->isAIJ    = PETSC_TRUE;
2600bf0cc555SLisandro Dalcin   mumps->Destroy  = B->ops->destroy;
2601450b117fSShri Abhyankar   B->ops->destroy = MatDestroy_MUMPS;
2602450b117fSShri Abhyankar   B->spptr        = (void*)mumps;
26032205254eSKarl Rupp 
2604f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
2605746480a1SHong Zhang 
2606450b117fSShri Abhyankar   *F = B;
2607450b117fSShri Abhyankar   PetscFunctionReturn(0);
2608450b117fSShri Abhyankar }
260942c9c57cSBarry Smith 
261042c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mumps(Mat,MatFactorType,Mat*);
261142c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_baij_mumps(Mat,MatFactorType,Mat*);
261242c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat,MatFactorType,Mat*);
261342c9c57cSBarry Smith 
261442c9c57cSBarry Smith #undef __FUNCT__
261542c9c57cSBarry Smith #define __FUNCT__ "MatSolverPackageRegister_MUMPS"
261629b38603SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_MUMPS(void)
261742c9c57cSBarry Smith {
261842c9c57cSBarry Smith   PetscErrorCode ierr;
261942c9c57cSBarry Smith 
262042c9c57cSBarry Smith   PetscFunctionBegin;
262142c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr);
262242c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr);
262342c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr);
262442c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr);
262542c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPISBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr);
262642c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr);
262742c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr);
262842c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr);
262942c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr);
263042c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr);
263142c9c57cSBarry Smith   PetscFunctionReturn(0);
263242c9c57cSBarry Smith }
263342c9c57cSBarry Smith 
2634