xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision 78a7176c4e3df89addbf1445c427b17e228bf17b)
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;
982205254eSKarl Rupp 
99bf0cc555SLisandro Dalcin   PetscErrorCode (*Destroy)(Mat);
100bccb9932SShri Abhyankar   PetscErrorCode (*ConvertToTriples)(Mat, int, MatReuse, int*, int**, int**, PetscScalar**);
101f0c56d0fSKris Buschelman } Mat_MUMPS;
102f0c56d0fSKris Buschelman 
10309573ac7SBarry Smith extern PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);
104b24902e0SBarry Smith 
10559ac8732SStefano Zampini #undef __FUNCT__
10659ac8732SStefano Zampini #define __FUNCT__ "MatMumpsResetSchur_Private"
10759ac8732SStefano Zampini static PetscErrorCode MatMumpsResetSchur_Private(Mat_MUMPS* mumps)
108b5fa320bSStefano Zampini {
109b5fa320bSStefano Zampini   PetscErrorCode ierr;
110b5fa320bSStefano Zampini 
111b5fa320bSStefano Zampini   PetscFunctionBegin;
11259ac8732SStefano Zampini   ierr = PetscFree2(mumps->id.listvar_schur,mumps->id.schur);CHKERRQ(ierr);
11359ac8732SStefano Zampini   ierr = PetscFree(mumps->id.redrhs);CHKERRQ(ierr);
11459ac8732SStefano Zampini   ierr = PetscFree(mumps->schur_sol);CHKERRQ(ierr);
11559ac8732SStefano Zampini   ierr = PetscFree(mumps->schur_pivots);CHKERRQ(ierr);
11659ac8732SStefano Zampini   ierr = PetscFree(mumps->schur_work);CHKERRQ(ierr);
11759ac8732SStefano Zampini   mumps->id.size_schur = 0;
11859ac8732SStefano Zampini   mumps->id.ICNTL(19) = 0;
11959ac8732SStefano Zampini   PetscFunctionReturn(0);
12059ac8732SStefano Zampini }
12159ac8732SStefano Zampini 
12259ac8732SStefano Zampini #undef __FUNCT__
12359ac8732SStefano Zampini #define __FUNCT__ "MatMumpsFactorSchur_Private"
12459ac8732SStefano Zampini static PetscErrorCode MatMumpsFactorSchur_Private(Mat_MUMPS* mumps)
12559ac8732SStefano Zampini {
12659ac8732SStefano Zampini   PetscBLASInt   B_N,B_ierr,B_slda;
12759ac8732SStefano Zampini   PetscErrorCode ierr;
12859ac8732SStefano Zampini 
12959ac8732SStefano Zampini   PetscFunctionBegin;
13059ac8732SStefano Zampini   if (mumps->schur_factored) {
13159ac8732SStefano Zampini     PetscFunctionReturn(0);
13259ac8732SStefano Zampini   }
13359ac8732SStefano Zampini   ierr = PetscBLASIntCast(mumps->id.size_schur,&B_N);CHKERRQ(ierr);
13459ac8732SStefano Zampini   ierr = PetscBLASIntCast(mumps->id.schur_lld,&B_slda);CHKERRQ(ierr);
13559ac8732SStefano Zampini   if (!mumps->sym) { /* MUMPS always return a full Schur matrix */
13659ac8732SStefano Zampini     if (!mumps->schur_pivots) {
13759ac8732SStefano Zampini       ierr = PetscMalloc1(B_N,&mumps->schur_pivots);CHKERRQ(ierr);
13859ac8732SStefano Zampini     }
13959ac8732SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
14059ac8732SStefano Zampini     PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,(PetscScalar*)mumps->id.schur,&B_slda,mumps->schur_pivots,&B_ierr));
14159ac8732SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
14259ac8732SStefano Zampini     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
14359ac8732SStefano Zampini   } else { /* either full or lower-triangular (not packed) */
14459ac8732SStefano Zampini     char ord[2];
14559ac8732SStefano Zampini     if (mumps->id.ICNTL(19) == 2 || mumps->id.ICNTL(19) == 3) { /* lower triangular stored by columns or full matrix */
14659ac8732SStefano Zampini       sprintf(ord,"L");
14759ac8732SStefano Zampini     } else { /* ICNTL(19) == 1 lower triangular stored by rows */
14859ac8732SStefano Zampini       sprintf(ord,"U");
14959ac8732SStefano Zampini     }
15059ac8732SStefano Zampini     if (mumps->id.sym == 2) {
15159ac8732SStefano Zampini       if (!mumps->schur_pivots) {
15259ac8732SStefano Zampini         PetscScalar  lwork;
15359ac8732SStefano Zampini 
15459ac8732SStefano Zampini         ierr = PetscMalloc1(B_N,&mumps->schur_pivots);CHKERRQ(ierr);
15559ac8732SStefano Zampini         mumps->schur_B_lwork=-1;
15659ac8732SStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
15759ac8732SStefano Zampini         PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_(ord,&B_N,(PetscScalar*)mumps->id.schur,&B_slda,mumps->schur_pivots,&lwork,&mumps->schur_B_lwork,&B_ierr));
15859ac8732SStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
15959ac8732SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYTRF Lapack routine %d",(int)B_ierr);
16059ac8732SStefano Zampini         ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&mumps->schur_B_lwork);CHKERRQ(ierr);
16159ac8732SStefano Zampini         ierr = PetscMalloc1(mumps->schur_B_lwork,&mumps->schur_work);CHKERRQ(ierr);
16259ac8732SStefano Zampini       }
16359ac8732SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
16459ac8732SStefano 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));
16559ac8732SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
16659ac8732SStefano Zampini       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
16759ac8732SStefano Zampini     } else {
16859ac8732SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
16959ac8732SStefano Zampini       PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_(ord,&B_N,(PetscScalar*)mumps->id.schur,&B_slda,&B_ierr));
17059ac8732SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
17159ac8732SStefano Zampini       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
17259ac8732SStefano Zampini     }
17359ac8732SStefano Zampini   }
17459ac8732SStefano Zampini   mumps->schur_factored = PETSC_TRUE;
17559ac8732SStefano Zampini   PetscFunctionReturn(0);
17659ac8732SStefano Zampini }
17759ac8732SStefano Zampini 
17859ac8732SStefano Zampini #undef __FUNCT__
17959ac8732SStefano Zampini #define __FUNCT__ "MatMumpsInvertSchur_Private"
18059ac8732SStefano Zampini static PetscErrorCode MatMumpsInvertSchur_Private(Mat_MUMPS* mumps)
18159ac8732SStefano Zampini {
18259ac8732SStefano Zampini   PetscBLASInt   B_N,B_ierr,B_slda;
18359ac8732SStefano Zampini   PetscErrorCode ierr;
18459ac8732SStefano Zampini 
18559ac8732SStefano Zampini   PetscFunctionBegin;
18659ac8732SStefano Zampini   ierr = MatMumpsFactorSchur_Private(mumps);CHKERRQ(ierr);
18759ac8732SStefano Zampini   ierr = PetscBLASIntCast(mumps->id.size_schur,&B_N);CHKERRQ(ierr);
18859ac8732SStefano Zampini   ierr = PetscBLASIntCast(mumps->id.schur_lld,&B_slda);CHKERRQ(ierr);
18959ac8732SStefano Zampini   if (!mumps->sym) { /* MUMPS always return a full Schur matrix */
19059ac8732SStefano Zampini     if (!mumps->schur_work) {
19159ac8732SStefano Zampini       PetscScalar lwork;
19259ac8732SStefano Zampini 
19359ac8732SStefano Zampini       mumps->schur_B_lwork = -1;
19459ac8732SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
19559ac8732SStefano Zampini       PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,(PetscScalar*)mumps->id.schur,&B_N,mumps->schur_pivots,&lwork,&mumps->schur_B_lwork,&B_ierr));
19659ac8732SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
19759ac8732SStefano Zampini       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
19859ac8732SStefano Zampini       ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&mumps->schur_B_lwork);CHKERRQ(ierr);
19959ac8732SStefano Zampini       ierr = PetscMalloc1(mumps->schur_B_lwork,&mumps->schur_work);CHKERRQ(ierr);
20059ac8732SStefano Zampini     }
20159ac8732SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
20259ac8732SStefano Zampini     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,(PetscScalar*)mumps->id.schur,&B_N,mumps->schur_pivots,mumps->schur_work,&mumps->schur_B_lwork,&B_ierr));
20359ac8732SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
20459ac8732SStefano Zampini     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
20559ac8732SStefano Zampini   } else { /* either full or lower-triangular (not packed) */
20659ac8732SStefano Zampini     char ord[2];
20759ac8732SStefano Zampini     if (mumps->id.ICNTL(19) == 2 || mumps->id.ICNTL(19) == 3) { /* lower triangular stored by columns or full matrix */
20859ac8732SStefano Zampini       sprintf(ord,"L");
20959ac8732SStefano Zampini     } else { /* ICNTL(19) == 1 lower triangular stored by rows */
21059ac8732SStefano Zampini       sprintf(ord,"U");
21159ac8732SStefano Zampini     }
21259ac8732SStefano Zampini     if (mumps->id.sym == 2) {
21359ac8732SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
21459ac8732SStefano Zampini       PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_(ord,&B_N,(PetscScalar*)mumps->id.schur,&B_N,mumps->schur_pivots,mumps->schur_work,&B_ierr));
21559ac8732SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
21659ac8732SStefano Zampini       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
21759ac8732SStefano Zampini     } else {
21859ac8732SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
21959ac8732SStefano Zampini       PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_(ord,&B_N,(PetscScalar*)mumps->id.schur,&B_N,&B_ierr));
22059ac8732SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
22159ac8732SStefano Zampini       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
22259ac8732SStefano Zampini     }
22359ac8732SStefano Zampini   }
22459ac8732SStefano Zampini   mumps->schur_inverted = PETSC_TRUE;
22559ac8732SStefano Zampini   PetscFunctionReturn(0);
22659ac8732SStefano Zampini }
22759ac8732SStefano Zampini 
22859ac8732SStefano Zampini #undef __FUNCT__
22959ac8732SStefano Zampini #define __FUNCT__ "MatMumpsSolveSchur_Private"
230e807eca7SStefano Zampini static PetscErrorCode MatMumpsSolveSchur_Private(Mat_MUMPS* mumps, PetscBool sol_in_redrhs)
23159ac8732SStefano Zampini {
23259ac8732SStefano Zampini   PetscBLASInt   B_N,B_Nrhs,B_ierr,B_slda,B_rlda;
23359ac8732SStefano Zampini   PetscScalar    one=1.,zero=0.;
23459ac8732SStefano Zampini   PetscErrorCode ierr;
23559ac8732SStefano Zampini 
23659ac8732SStefano Zampini   PetscFunctionBegin;
23759ac8732SStefano Zampini   ierr = MatMumpsFactorSchur_Private(mumps);CHKERRQ(ierr);
238b5fa320bSStefano Zampini   ierr = PetscBLASIntCast(mumps->id.size_schur,&B_N);CHKERRQ(ierr);
239b5fa320bSStefano Zampini   ierr = PetscBLASIntCast(mumps->id.schur_lld,&B_slda);CHKERRQ(ierr);
240b5fa320bSStefano Zampini   ierr = PetscBLASIntCast(mumps->id.nrhs,&B_Nrhs);CHKERRQ(ierr);
241b5fa320bSStefano Zampini   ierr = PetscBLASIntCast(mumps->id.lredrhs,&B_rlda);CHKERRQ(ierr);
24259ac8732SStefano Zampini   if (mumps->schur_inverted) {
24359ac8732SStefano Zampini     PetscInt sizesol = B_Nrhs*B_N;
24459ac8732SStefano Zampini     if (!mumps->schur_sol || sizesol > mumps->schur_sizesol) {
24559ac8732SStefano Zampini       ierr = PetscFree(mumps->schur_sol);CHKERRQ(ierr);
24659ac8732SStefano Zampini       ierr = PetscMalloc1(sizesol,&mumps->schur_sol);CHKERRQ(ierr);
24759ac8732SStefano Zampini       mumps->schur_sizesol = sizesol;
248b5fa320bSStefano Zampini     }
24959ac8732SStefano Zampini     if (!mumps->sym) {
25059ac8732SStefano Zampini       char type[2];
251b5fa320bSStefano Zampini       if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
25259ac8732SStefano Zampini         if (!mumps->id.ICNTL(9)) { /* transpose solve */
25359ac8732SStefano Zampini           sprintf(type,"N");
254b5fa320bSStefano Zampini         } else {
25559ac8732SStefano Zampini           sprintf(type,"T");
256b5fa320bSStefano Zampini         }
25759ac8732SStefano Zampini       } else { /* stored by columns */
25859ac8732SStefano Zampini         if (!mumps->id.ICNTL(9)) { /* transpose solve */
25959ac8732SStefano Zampini           sprintf(type,"T");
26059ac8732SStefano Zampini         } else {
26159ac8732SStefano Zampini           sprintf(type,"N");
26259ac8732SStefano Zampini         }
26359ac8732SStefano Zampini       }
26459ac8732SStefano 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));
26559ac8732SStefano Zampini     } else {
26659ac8732SStefano Zampini       char ord[2];
26759ac8732SStefano Zampini       if (mumps->id.ICNTL(19) == 2 || mumps->id.ICNTL(19) == 3) { /* lower triangular stored by columns or full matrix */
26859ac8732SStefano Zampini         sprintf(ord,"L");
26959ac8732SStefano Zampini       } else { /* ICNTL(19) == 1 lower triangular stored by rows */
27059ac8732SStefano Zampini         sprintf(ord,"U");
27159ac8732SStefano Zampini       }
27259ac8732SStefano 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));
27359ac8732SStefano Zampini     }
274e807eca7SStefano Zampini     if (sol_in_redrhs) {
27559ac8732SStefano Zampini       ierr = PetscMemcpy(mumps->id.redrhs,mumps->schur_sol,sizesol*sizeof(PetscScalar));CHKERRQ(ierr);
276e807eca7SStefano Zampini     }
277a12f35bfSStefano Zampini   } else { /* Schur complement has not been inverted */
278a12f35bfSStefano Zampini     MumpsScalar *orhs=NULL;
279a12f35bfSStefano Zampini 
280a12f35bfSStefano Zampini     if (!sol_in_redrhs) {
281a12f35bfSStefano Zampini       PetscInt sizesol = B_Nrhs*B_N;
282a12f35bfSStefano Zampini       if (!mumps->schur_sol || sizesol > mumps->schur_sizesol) {
283a12f35bfSStefano Zampini         ierr = PetscFree(mumps->schur_sol);CHKERRQ(ierr);
284a12f35bfSStefano Zampini         ierr = PetscMalloc1(sizesol,&mumps->schur_sol);CHKERRQ(ierr);
285a12f35bfSStefano Zampini         mumps->schur_sizesol = sizesol;
286a12f35bfSStefano Zampini       }
287a12f35bfSStefano Zampini       orhs = mumps->id.redrhs;
288a12f35bfSStefano Zampini       ierr = PetscMemcpy(mumps->schur_sol,mumps->id.redrhs,sizesol*sizeof(PetscScalar));CHKERRQ(ierr);
289a12f35bfSStefano Zampini       mumps->id.redrhs = (MumpsScalar*)mumps->schur_sol;
290a12f35bfSStefano Zampini     }
29159ac8732SStefano Zampini     if (!mumps->sym) { /* MUMPS always return a full Schur matrix */
29259ac8732SStefano Zampini       char type[2];
29359ac8732SStefano Zampini       if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
29459ac8732SStefano Zampini         if (!mumps->id.ICNTL(9)) { /* transpose solve */
29559ac8732SStefano Zampini           sprintf(type,"N");
29659ac8732SStefano Zampini         } else {
29759ac8732SStefano Zampini           sprintf(type,"T");
29859ac8732SStefano Zampini         }
29959ac8732SStefano Zampini       } else { /* stored by columns */
30059ac8732SStefano Zampini         if (!mumps->id.ICNTL(9)) { /* transpose solve */
30159ac8732SStefano Zampini           sprintf(type,"T");
30259ac8732SStefano Zampini         } else {
30359ac8732SStefano Zampini           sprintf(type,"N");
30459ac8732SStefano Zampini         }
30559ac8732SStefano Zampini       }
30659ac8732SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
30759ac8732SStefano 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));
308b5fa320bSStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
309b5fa320bSStefano Zampini       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRS Lapack routine %d",(int)B_ierr);
310b5fa320bSStefano Zampini     } else { /* either full or lower-triangular (not packed) */
311b5fa320bSStefano Zampini       char ord[2];
312b5fa320bSStefano Zampini       if (mumps->id.ICNTL(19) == 2 || mumps->id.ICNTL(19) == 3) { /* lower triangular stored by columns or full matrix */
313b5fa320bSStefano Zampini         sprintf(ord,"L");
314b5fa320bSStefano Zampini       } else { /* ICNTL(19) == 1 lower triangular stored by rows */
315b5fa320bSStefano Zampini         sprintf(ord,"U");
316b5fa320bSStefano Zampini       }
317b5fa320bSStefano Zampini       if (mumps->id.sym == 2) {
318b5fa320bSStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
31959ac8732SStefano 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));
320b5fa320bSStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
321b5fa320bSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRS Lapack routine %d",(int)B_ierr);
322b5fa320bSStefano Zampini       } else {
323b5fa320bSStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
32459ac8732SStefano Zampini         PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_(ord,&B_N,&B_Nrhs,(PetscScalar*)mumps->id.schur,&B_slda,(PetscScalar*)mumps->id.redrhs,&B_rlda,&B_ierr));
325b5fa320bSStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
326b5fa320bSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRS Lapack routine %d",(int)B_ierr);
327b5fa320bSStefano Zampini       }
328b5fa320bSStefano Zampini     }
329e807eca7SStefano Zampini     if (!sol_in_redrhs) {
330a12f35bfSStefano Zampini       mumps->id.redrhs = orhs;
331e807eca7SStefano Zampini     }
33259ac8732SStefano Zampini   }
333b5fa320bSStefano Zampini   PetscFunctionReturn(0);
334b5fa320bSStefano Zampini }
335b5fa320bSStefano Zampini 
33659ac8732SStefano Zampini #undef __FUNCT__
33759ac8732SStefano Zampini #define __FUNCT__ "MatMumpsHandleSchur_Private"
338b8f61ee1SStefano Zampini static PetscErrorCode MatMumpsHandleSchur_Private(Mat_MUMPS* mumps, PetscBool expansion)
339b5fa320bSStefano Zampini {
340b5fa320bSStefano Zampini   PetscErrorCode ierr;
341b5fa320bSStefano Zampini 
342b5fa320bSStefano Zampini   PetscFunctionBegin;
343b5fa320bSStefano Zampini   if (!mumps->id.ICNTL(19)) { /* do nothing when Schur complement has not been computed */
344b5fa320bSStefano Zampini     PetscFunctionReturn(0);
345b5fa320bSStefano Zampini   }
346b8f61ee1SStefano Zampini   if (!expansion) { /* prepare for the condensation step */
347b5fa320bSStefano Zampini     /* check if schur complement has been computed
348e807eca7SStefano Zampini        We set by default ICNTL(26) == -1 when Schur indices have been provided by the user.
349b5fa320bSStefano Zampini        According to MUMPS (5.0.0) manual, any value should be harmful during the factorization phase
350b5fa320bSStefano Zampini        Unless the user provides a valid value for ICNTL(26), MatSolve and MatMatSolve routines solve the full system.
351b8f61ee1SStefano Zampini        This requires an extra call to PetscMUMPS_c and the computation of the factors for S */
352b5fa320bSStefano Zampini     if (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2) {
353b5fa320bSStefano Zampini       PetscInt sizeredrhs = mumps->id.nrhs*mumps->id.size_schur;
354b5fa320bSStefano Zampini       /* allocate MUMPS internal array to store reduced right-hand sides */
355b5fa320bSStefano Zampini       if (!mumps->id.redrhs || sizeredrhs > mumps->sizeredrhs) {
356b5fa320bSStefano Zampini         ierr = PetscFree(mumps->id.redrhs);CHKERRQ(ierr);
357b5fa320bSStefano Zampini         mumps->id.lredrhs = mumps->id.size_schur;
358b5fa320bSStefano Zampini         ierr = PetscMalloc1(mumps->id.nrhs*mumps->id.lredrhs,&mumps->id.redrhs);CHKERRQ(ierr);
359b5fa320bSStefano Zampini         mumps->sizeredrhs = mumps->id.nrhs*mumps->id.lredrhs;
360b5fa320bSStefano Zampini       }
361b5fa320bSStefano Zampini       mumps->id.ICNTL(26) = 1; /* condensation phase */
362b5fa320bSStefano Zampini     }
363b5fa320bSStefano Zampini   } else { /* prepare for the expansion step */
364b8f61ee1SStefano Zampini     /* solve Schur complement (this has to be done by the MUMPS user, so basically us) */
365e807eca7SStefano Zampini     ierr = MatMumpsSolveSchur_Private(mumps,PETSC_TRUE);CHKERRQ(ierr);
366b5fa320bSStefano Zampini     mumps->id.ICNTL(26) = 2; /* expansion phase */
367b5fa320bSStefano Zampini     PetscMUMPS_c(&mumps->id);
368b5fa320bSStefano 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));
369b5fa320bSStefano Zampini     /* restore defaults */
370b5fa320bSStefano Zampini     mumps->id.ICNTL(26) = -1;
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);
857bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetIcntl_C",NULL);CHKERRQ(ierr);
858bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetIcntl_C",NULL);CHKERRQ(ierr);
859bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetCntl_C",NULL);CHKERRQ(ierr);
860bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetCntl_C",NULL);CHKERRQ(ierr);
861ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfo_C",NULL);CHKERRQ(ierr);
862ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfog_C",NULL);CHKERRQ(ierr);
863ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfo_C",NULL);CHKERRQ(ierr);
864ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfog_C",NULL);CHKERRQ(ierr);
865397b6df1SKris Buschelman   PetscFunctionReturn(0);
866397b6df1SKris Buschelman }
867397b6df1SKris Buschelman 
868397b6df1SKris Buschelman #undef __FUNCT__
869f6c57405SHong Zhang #define __FUNCT__ "MatSolve_MUMPS"
870b24902e0SBarry Smith PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x)
871b24902e0SBarry Smith {
872a5e57a09SHong Zhang   Mat_MUMPS        *mumps=(Mat_MUMPS*)A->spptr;
873d54de34fSKris Buschelman   PetscScalar      *array;
87467877ebaSShri Abhyankar   Vec              b_seq;
875329ec9b3SHong Zhang   IS               is_iden,is_petsc;
876dfbe8321SBarry Smith   PetscErrorCode   ierr;
877329ec9b3SHong Zhang   PetscInt         i;
878883f2eb9SBarry Smith   static PetscBool cite1 = PETSC_FALSE,cite2 = PETSC_FALSE;
879397b6df1SKris Buschelman 
880397b6df1SKris Buschelman   PetscFunctionBegin;
881883f2eb9SBarry 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);
882883f2eb9SBarry 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);
883a5e57a09SHong Zhang   mumps->id.nrhs = 1;
884a5e57a09SHong Zhang   b_seq          = mumps->b_seq;
885a5e57a09SHong Zhang   if (mumps->size > 1) {
886329ec9b3SHong Zhang     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
887a5e57a09SHong Zhang     ierr = VecScatterBegin(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
888a5e57a09SHong Zhang     ierr = VecScatterEnd(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
889a5e57a09SHong Zhang     if (!mumps->myid) {ierr = VecGetArray(b_seq,&array);CHKERRQ(ierr);}
890397b6df1SKris Buschelman   } else {  /* size == 1 */
891397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
892397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
893397b6df1SKris Buschelman   }
894a5e57a09SHong Zhang   if (!mumps->myid) { /* define rhs on the host */
895a5e57a09SHong Zhang     mumps->id.nrhs = 1;
896940cd9d6SSatish Balay     mumps->id.rhs = (MumpsScalar*)array;
897397b6df1SKris Buschelman   }
898397b6df1SKris Buschelman 
899b5fa320bSStefano Zampini   /* handle condensation step of Schur complement (if any) */
900b8f61ee1SStefano Zampini   ierr = MatMumpsHandleSchur_Private(mumps,PETSC_FALSE);CHKERRQ(ierr);
901b5fa320bSStefano Zampini 
902397b6df1SKris Buschelman   /* solve phase */
903329ec9b3SHong Zhang   /*-------------*/
904a5e57a09SHong Zhang   mumps->id.job = JOB_SOLVE;
905a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
906a5e57a09SHong 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));
907397b6df1SKris Buschelman 
908b5fa320bSStefano Zampini   /* handle expansion step of Schur complement (if any) */
909b8f61ee1SStefano Zampini   ierr = MatMumpsHandleSchur_Private(mumps,PETSC_TRUE);CHKERRQ(ierr);
910b5fa320bSStefano Zampini 
911a5e57a09SHong Zhang   if (mumps->size > 1) { /* convert mumps distributed solution to petsc mpi x */
912a5e57a09SHong Zhang     if (mumps->scat_sol && mumps->ICNTL9_pre != mumps->id.ICNTL(9)) {
913a5e57a09SHong Zhang       /* when id.ICNTL(9) changes, the contents of lsol_loc may change (not its size, lsol_loc), recreates scat_sol */
914a5e57a09SHong Zhang       ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
915397b6df1SKris Buschelman     }
916a5e57a09SHong Zhang     if (!mumps->scat_sol) { /* create scatter scat_sol */
917a5e57a09SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,mumps->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */
918a5e57a09SHong Zhang       for (i=0; i<mumps->id.lsol_loc; i++) {
919a5e57a09SHong Zhang         mumps->id.isol_loc[i] -= 1; /* change Fortran style to C style */
920a5e57a09SHong Zhang       }
921a5e57a09SHong Zhang       ierr = ISCreateGeneral(PETSC_COMM_SELF,mumps->id.lsol_loc,mumps->id.isol_loc,PETSC_COPY_VALUES,&is_petsc);CHKERRQ(ierr);  /* to */
922a5e57a09SHong Zhang       ierr = VecScatterCreate(mumps->x_seq,is_iden,x,is_petsc,&mumps->scat_sol);CHKERRQ(ierr);
9236bf464f9SBarry Smith       ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
9246bf464f9SBarry Smith       ierr = ISDestroy(&is_petsc);CHKERRQ(ierr);
9252205254eSKarl Rupp 
926a5e57a09SHong Zhang       mumps->ICNTL9_pre = mumps->id.ICNTL(9); /* save current value of id.ICNTL(9) */
927397b6df1SKris Buschelman     }
928a5e57a09SHong Zhang 
929a5e57a09SHong Zhang     ierr = VecScatterBegin(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
930a5e57a09SHong Zhang     ierr = VecScatterEnd(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
931329ec9b3SHong Zhang   }
932397b6df1SKris Buschelman   PetscFunctionReturn(0);
933397b6df1SKris Buschelman }
934397b6df1SKris Buschelman 
93551d5961aSHong Zhang #undef __FUNCT__
93651d5961aSHong Zhang #define __FUNCT__ "MatSolveTranspose_MUMPS"
93751d5961aSHong Zhang PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x)
93851d5961aSHong Zhang {
939a5e57a09SHong Zhang   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
94051d5961aSHong Zhang   PetscErrorCode ierr;
94151d5961aSHong Zhang 
94251d5961aSHong Zhang   PetscFunctionBegin;
943a5e57a09SHong Zhang   mumps->id.ICNTL(9) = 0;
9440ad0caddSJed Brown   ierr = MatSolve_MUMPS(A,b,x);CHKERRQ(ierr);
945a5e57a09SHong Zhang   mumps->id.ICNTL(9) = 1;
94651d5961aSHong Zhang   PetscFunctionReturn(0);
94751d5961aSHong Zhang }
94851d5961aSHong Zhang 
949e0b74bf9SHong Zhang #undef __FUNCT__
950e0b74bf9SHong Zhang #define __FUNCT__ "MatMatSolve_MUMPS"
951e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X)
952e0b74bf9SHong Zhang {
953bda8bf91SBarry Smith   PetscErrorCode ierr;
954bda8bf91SBarry Smith   PetscBool      flg;
9554e34a73bSHong Zhang   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
956334c5f61SHong Zhang   PetscInt       i,nrhs,M;
9572cd7d884SHong Zhang   PetscScalar    *array,*bray;
958bda8bf91SBarry Smith 
959e0b74bf9SHong Zhang   PetscFunctionBegin;
9600298fd71SBarry Smith   ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
961801fbe65SHong Zhang   if (!flg) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
9620298fd71SBarry Smith   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
963801fbe65SHong Zhang   if (!flg) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
964801fbe65SHong 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");
9654e34a73bSHong Zhang 
9662cd7d884SHong Zhang   ierr = MatGetSize(B,&M,&nrhs);CHKERRQ(ierr);
967334c5f61SHong Zhang   mumps->id.nrhs = nrhs;
968334c5f61SHong Zhang   mumps->id.lrhs = M;
9694e34a73bSHong Zhang 
9702cd7d884SHong Zhang   if (mumps->size == 1) {
9712cd7d884SHong Zhang     /* copy B to X */
9722cd7d884SHong Zhang     ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr);
9732cd7d884SHong Zhang     ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr);
9746444a565SStefano Zampini     ierr = PetscMemcpy(array,bray,M*nrhs*sizeof(PetscScalar));CHKERRQ(ierr);
9752cd7d884SHong Zhang     ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr);
976940cd9d6SSatish Balay     mumps->id.rhs = (MumpsScalar*)array;
977b5fa320bSStefano Zampini     /* handle condensation step of Schur complement (if any) */
978b8f61ee1SStefano Zampini     ierr = MatMumpsHandleSchur_Private(mumps,PETSC_FALSE);CHKERRQ(ierr);
979801fbe65SHong Zhang 
9802cd7d884SHong Zhang     /* solve phase */
9812cd7d884SHong Zhang     /*-------------*/
9822cd7d884SHong Zhang     mumps->id.job = JOB_SOLVE;
9832cd7d884SHong Zhang     PetscMUMPS_c(&mumps->id);
9842cd7d884SHong 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));
985b5fa320bSStefano Zampini 
986b5fa320bSStefano Zampini     /* handle expansion step of Schur complement (if any) */
987b8f61ee1SStefano Zampini     ierr = MatMumpsHandleSchur_Private(mumps,PETSC_TRUE);CHKERRQ(ierr);
9882cd7d884SHong Zhang     ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr);
989334c5f61SHong Zhang   } else {  /*--------- parallel case --------*/
99071aed81dSHong Zhang     PetscInt       lsol_loc,nlsol_loc,*isol_loc,*idx,*iidx,*idxx,*isol_loc_save;
9911070efccSSatish Balay     MumpsScalar    *sol_loc,*sol_loc_save;
992801fbe65SHong Zhang     IS             is_to,is_from;
993334c5f61SHong Zhang     PetscInt       k,proc,j,m;
994801fbe65SHong Zhang     const PetscInt *rstart;
995334c5f61SHong Zhang     Vec            v_mpi,b_seq,x_seq;
996334c5f61SHong Zhang     VecScatter     scat_rhs,scat_sol;
997801fbe65SHong Zhang 
998801fbe65SHong Zhang     /* create x_seq to hold local solution */
99971aed81dSHong Zhang     isol_loc_save = mumps->id.isol_loc; /* save it for MatSovle() */
100071aed81dSHong Zhang     sol_loc_save  = mumps->id.sol_loc;
1001801fbe65SHong Zhang 
100271aed81dSHong Zhang     lsol_loc  = mumps->id.INFO(23);
100371aed81dSHong Zhang     nlsol_loc = nrhs*lsol_loc;     /* length of sol_loc */
100471aed81dSHong Zhang     ierr = PetscMalloc2(nlsol_loc,&sol_loc,nlsol_loc,&isol_loc);CHKERRQ(ierr);
1005940cd9d6SSatish Balay     mumps->id.sol_loc = (MumpsScalar*)sol_loc;
1006801fbe65SHong Zhang     mumps->id.isol_loc = isol_loc;
1007801fbe65SHong Zhang 
10081070efccSSatish Balay     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,nlsol_loc,(PetscScalar*)sol_loc,&x_seq);CHKERRQ(ierr);
10092cd7d884SHong Zhang 
101074f0fcc7SHong Zhang     /* copy rhs matrix B into vector v_mpi */
1011334c5f61SHong Zhang     ierr = MatGetLocalSize(B,&m,NULL);CHKERRQ(ierr);
1012801fbe65SHong Zhang     ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr);
101374f0fcc7SHong Zhang     ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)B),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi);CHKERRQ(ierr);
1014801fbe65SHong Zhang     ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr);
1015801fbe65SHong Zhang 
1016334c5f61SHong Zhang     /* scatter v_mpi to b_seq because MUMPS only supports centralized rhs */
101774f0fcc7SHong Zhang     /* idx: maps from k-th index of v_mpi to (i,j)-th global entry of B;
1018801fbe65SHong Zhang       iidx: inverse of idx, will be used by scattering xx_seq -> X       */
1019801fbe65SHong Zhang     ierr = PetscMalloc2(nrhs*M,&idx,nrhs*M,&iidx);CHKERRQ(ierr);
1020801fbe65SHong Zhang     ierr = MatGetOwnershipRanges(B,&rstart);CHKERRQ(ierr);
1021801fbe65SHong Zhang     k = 0;
1022801fbe65SHong Zhang     for (proc=0; proc<mumps->size; proc++){
1023801fbe65SHong Zhang       for (j=0; j<nrhs; j++){
1024801fbe65SHong Zhang         for (i=rstart[proc]; i<rstart[proc+1]; i++){
1025801fbe65SHong Zhang           iidx[j*M + i] = k;
1026801fbe65SHong Zhang           idx[k++]      = j*M + i;
1027801fbe65SHong Zhang         }
1028801fbe65SHong Zhang       }
10292cd7d884SHong Zhang     }
10302cd7d884SHong Zhang 
1031801fbe65SHong Zhang     if (!mumps->myid) {
1032334c5f61SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,nrhs*M,&b_seq);CHKERRQ(ierr);
1033801fbe65SHong Zhang       ierr = ISCreateGeneral(PETSC_COMM_SELF,nrhs*M,idx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr);
1034801fbe65SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,nrhs*M,0,1,&is_from);CHKERRQ(ierr);
1035801fbe65SHong Zhang     } else {
1036334c5f61SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&b_seq);CHKERRQ(ierr);
1037801fbe65SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_to);CHKERRQ(ierr);
1038801fbe65SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_from);CHKERRQ(ierr);
1039801fbe65SHong Zhang     }
1040334c5f61SHong Zhang     ierr = VecScatterCreate(v_mpi,is_from,b_seq,is_to,&scat_rhs);CHKERRQ(ierr);
1041334c5f61SHong Zhang     ierr = VecScatterBegin(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1042801fbe65SHong Zhang     ierr = ISDestroy(&is_to);CHKERRQ(ierr);
1043801fbe65SHong Zhang     ierr = ISDestroy(&is_from);CHKERRQ(ierr);
1044334c5f61SHong Zhang     ierr = VecScatterEnd(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1045801fbe65SHong Zhang 
1046801fbe65SHong Zhang     if (!mumps->myid) { /* define rhs on the host */
1047334c5f61SHong Zhang       ierr = VecGetArray(b_seq,&bray);CHKERRQ(ierr);
1048940cd9d6SSatish Balay       mumps->id.rhs = (MumpsScalar*)bray;
1049334c5f61SHong Zhang       ierr = VecRestoreArray(b_seq,&bray);CHKERRQ(ierr);
1050801fbe65SHong Zhang     }
1051801fbe65SHong Zhang 
1052801fbe65SHong Zhang     /* solve phase */
1053801fbe65SHong Zhang     /*-------------*/
1054801fbe65SHong Zhang     mumps->id.job = JOB_SOLVE;
1055801fbe65SHong Zhang     PetscMUMPS_c(&mumps->id);
1056801fbe65SHong 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));
1057801fbe65SHong Zhang 
1058334c5f61SHong Zhang     /* scatter mumps distributed solution to petsc vector v_mpi, which shares local arrays with solution matrix X */
105974f0fcc7SHong Zhang     ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr);
106074f0fcc7SHong Zhang     ierr = VecPlaceArray(v_mpi,array);CHKERRQ(ierr);
1061801fbe65SHong Zhang 
1062334c5f61SHong Zhang     /* create scatter scat_sol */
106371aed81dSHong Zhang     ierr = PetscMalloc1(nlsol_loc,&idxx);CHKERRQ(ierr);
106471aed81dSHong Zhang     ierr = ISCreateStride(PETSC_COMM_SELF,nlsol_loc,0,1,&is_from);CHKERRQ(ierr);
106571aed81dSHong Zhang     for (i=0; i<lsol_loc; i++) {
1066334c5f61SHong Zhang       isol_loc[i] -= 1; /* change Fortran style to C style */
1067334c5f61SHong Zhang       idxx[i] = iidx[isol_loc[i]];
1068801fbe65SHong Zhang       for (j=1; j<nrhs; j++){
1069334c5f61SHong Zhang         idxx[j*lsol_loc+i] = iidx[isol_loc[i]+j*M];
1070801fbe65SHong Zhang       }
1071801fbe65SHong Zhang     }
107271aed81dSHong Zhang     ierr = ISCreateGeneral(PETSC_COMM_SELF,nlsol_loc,idxx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr);
1073334c5f61SHong Zhang     ierr = VecScatterCreate(x_seq,is_from,v_mpi,is_to,&scat_sol);CHKERRQ(ierr);
1074334c5f61SHong Zhang     ierr = VecScatterBegin(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1075801fbe65SHong Zhang     ierr = ISDestroy(&is_from);CHKERRQ(ierr);
1076801fbe65SHong Zhang     ierr = ISDestroy(&is_to);CHKERRQ(ierr);
1077334c5f61SHong Zhang     ierr = VecScatterEnd(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1078801fbe65SHong Zhang     ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr);
107971aed81dSHong Zhang 
108071aed81dSHong Zhang     /* free spaces */
108171aed81dSHong Zhang     mumps->id.sol_loc = sol_loc_save;
108271aed81dSHong Zhang     mumps->id.isol_loc = isol_loc_save;
108371aed81dSHong Zhang 
108471aed81dSHong Zhang     ierr = PetscFree2(sol_loc,isol_loc);CHKERRQ(ierr);
1085801fbe65SHong Zhang     ierr = PetscFree2(idx,iidx);CHKERRQ(ierr);
1086801fbe65SHong Zhang     ierr = PetscFree(idxx);CHKERRQ(ierr);
108771aed81dSHong Zhang     ierr = VecDestroy(&x_seq);CHKERRQ(ierr);
108874f0fcc7SHong Zhang     ierr = VecDestroy(&v_mpi);CHKERRQ(ierr);
1089334c5f61SHong Zhang     ierr = VecDestroy(&b_seq);CHKERRQ(ierr);
1090334c5f61SHong Zhang     ierr = VecScatterDestroy(&scat_rhs);CHKERRQ(ierr);
1091334c5f61SHong Zhang     ierr = VecScatterDestroy(&scat_sol);CHKERRQ(ierr);
1092801fbe65SHong Zhang   }
1093e0b74bf9SHong Zhang   PetscFunctionReturn(0);
1094e0b74bf9SHong Zhang }
1095e0b74bf9SHong Zhang 
1096ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX)
1097a58c3f20SHong Zhang /*
1098a58c3f20SHong Zhang   input:
1099a58c3f20SHong Zhang    F:        numeric factor
1100a58c3f20SHong Zhang   output:
1101a58c3f20SHong Zhang    nneg:     total number of negative pivots
1102a58c3f20SHong Zhang    nzero:    0
1103a58c3f20SHong Zhang    npos:     (global dimension of F) - nneg
1104a58c3f20SHong Zhang */
1105a58c3f20SHong Zhang 
1106a58c3f20SHong Zhang #undef __FUNCT__
1107a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
1108dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
1109a58c3f20SHong Zhang {
1110a5e57a09SHong Zhang   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
1111dfbe8321SBarry Smith   PetscErrorCode ierr;
1112c1490034SHong Zhang   PetscMPIInt    size;
1113a58c3f20SHong Zhang 
1114a58c3f20SHong Zhang   PetscFunctionBegin;
1115ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)F),&size);CHKERRQ(ierr);
1116bcb30aebSHong 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 */
1117a5e57a09SHong 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));
1118ed85ac9fSHong Zhang 
1119710ac8efSHong Zhang   if (nneg) *nneg = mumps->id.INFOG(12);
1120ed85ac9fSHong Zhang   if (nzero || npos) {
1121ed85ac9fSHong 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");
1122710ac8efSHong Zhang     if (nzero) *nzero = mumps->id.INFOG(28);
1123710ac8efSHong Zhang     if (npos) *npos   = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28));
1124a58c3f20SHong Zhang   }
1125a58c3f20SHong Zhang   PetscFunctionReturn(0);
1126a58c3f20SHong Zhang }
1127ace3df97SHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */
1128a58c3f20SHong Zhang 
1129397b6df1SKris Buschelman #undef __FUNCT__
1130f6c57405SHong Zhang #define __FUNCT__ "MatFactorNumeric_MUMPS"
11310481f469SBarry Smith PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
1132af281ebdSHong Zhang {
1133a5e57a09SHong Zhang   Mat_MUMPS      *mumps =(Mat_MUMPS*)(F)->spptr;
11346849ba73SBarry Smith   PetscErrorCode ierr;
1135e09efc27SHong Zhang   Mat            F_diag;
1136ace3abfcSBarry Smith   PetscBool      isMPIAIJ;
1137397b6df1SKris Buschelman 
1138397b6df1SKris Buschelman   PetscFunctionBegin;
1139a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1140397b6df1SKris Buschelman 
1141397b6df1SKris Buschelman   /* numerical factorization phase */
1142329ec9b3SHong Zhang   /*-------------------------------*/
1143a5e57a09SHong Zhang   mumps->id.job = JOB_FACTNUMERIC;
11444e34a73bSHong Zhang   if (!mumps->id.ICNTL(18)) { /* A is centralized */
1145a5e57a09SHong Zhang     if (!mumps->myid) {
1146940cd9d6SSatish Balay       mumps->id.a = (MumpsScalar*)mumps->val;
1147397b6df1SKris Buschelman     }
1148397b6df1SKris Buschelman   } else {
1149940cd9d6SSatish Balay     mumps->id.a_loc = (MumpsScalar*)mumps->val;
1150397b6df1SKris Buschelman   }
1151a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
1152a5e57a09SHong Zhang   if (mumps->id.INFOG(1) < 0) {
1153151787a6SHong Zhang     if (mumps->id.INFO(1) == -13) {
1154151787a6SHong Zhang       if (mumps->id.INFO(2) < 0) {
1155151787a6SHong 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));
1156151787a6SHong Zhang       } else {
1157151787a6SHong 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));
1158151787a6SHong Zhang       }
1159151787a6SHong 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));
1160397b6df1SKris Buschelman   }
1161a5e57a09SHong 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));
1162397b6df1SKris Buschelman 
1163dcd589f8SShri Abhyankar   (F)->assembled        = PETSC_TRUE;
1164a5e57a09SHong Zhang   mumps->matstruc       = SAME_NONZERO_PATTERN;
1165b5fa320bSStefano Zampini   mumps->schur_factored = PETSC_FALSE;
116659ac8732SStefano Zampini   mumps->schur_inverted = PETSC_FALSE;
116767877ebaSShri Abhyankar 
1168066565c5SStefano Zampini   /* just to be sure that ICNTL(19) value returned by a call from MatMumpsGetIcntl is always consistent */
1169066565c5SStefano Zampini   if (!mumps->sym && mumps->id.ICNTL(19) && mumps->id.ICNTL(19) != 1) mumps->id.ICNTL(19) = 3;
1170066565c5SStefano Zampini 
1171a5e57a09SHong Zhang   if (mumps->size > 1) {
117267877ebaSShri Abhyankar     PetscInt    lsol_loc;
117367877ebaSShri Abhyankar     PetscScalar *sol_loc;
11742205254eSKarl Rupp 
1175c2093ab7SHong Zhang     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr);
1176c2093ab7SHong Zhang     if (isMPIAIJ) F_diag = ((Mat_MPIAIJ*)(F)->data)->A;
1177c2093ab7SHong Zhang     else F_diag = ((Mat_MPISBAIJ*)(F)->data)->A;
1178c2093ab7SHong Zhang     F_diag->assembled = PETSC_TRUE;
1179c2093ab7SHong Zhang 
1180c2093ab7SHong Zhang     /* distributed solution; Create x_seq=sol_loc for repeated use */
1181c2093ab7SHong Zhang     if (mumps->x_seq) {
1182c2093ab7SHong Zhang       ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
1183c2093ab7SHong Zhang       ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr);
1184c2093ab7SHong Zhang       ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr);
1185c2093ab7SHong Zhang     }
1186a5e57a09SHong Zhang     lsol_loc = mumps->id.INFO(23); /* length of sol_loc */
1187dcca6d9dSJed Brown     ierr = PetscMalloc2(lsol_loc,&sol_loc,lsol_loc,&mumps->id.isol_loc);CHKERRQ(ierr);
1188a5e57a09SHong Zhang     mumps->id.lsol_loc = lsol_loc;
1189940cd9d6SSatish Balay     mumps->id.sol_loc = (MumpsScalar*)sol_loc;
1190a5e57a09SHong Zhang     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq);CHKERRQ(ierr);
119167877ebaSShri Abhyankar   }
1192397b6df1SKris Buschelman   PetscFunctionReturn(0);
1193397b6df1SKris Buschelman }
1194397b6df1SKris Buschelman 
11959a2535b5SHong Zhang /* Sets MUMPS options from the options database */
1196dcd589f8SShri Abhyankar #undef __FUNCT__
11979a2535b5SHong Zhang #define __FUNCT__ "PetscSetMUMPSFromOptions"
11989a2535b5SHong Zhang PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A)
1199dcd589f8SShri Abhyankar {
12009a2535b5SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
1201dcd589f8SShri Abhyankar   PetscErrorCode ierr;
1202b34f08ffSHong Zhang   PetscInt       icntl,info[40],i,ninfo=40;
1203ace3abfcSBarry Smith   PetscBool      flg;
1204dcd589f8SShri Abhyankar 
1205dcd589f8SShri Abhyankar   PetscFunctionBegin;
1206ce94432eSBarry Smith   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
12079a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);CHKERRQ(ierr);
12089a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(1) = icntl;
12099a2535b5SHong 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);
12109a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(2) = icntl;
12119a2535b5SHong 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);
12129a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(3) = icntl;
1213dcd589f8SShri Abhyankar 
12149a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",mumps->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
12159a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(4) = icntl;
12169a2535b5SHong Zhang   if (mumps->id.ICNTL(4) || PetscLogPrintInfo) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */
12179a2535b5SHong Zhang 
1218d341cd04SHong 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);
12199a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(6) = icntl;
12209a2535b5SHong Zhang 
1221d341cd04SHong 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);
1222dcd589f8SShri Abhyankar   if (flg) {
12232205254eSKarl 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");
12242205254eSKarl Rupp     else mumps->id.ICNTL(7) = icntl;
1225dcd589f8SShri Abhyankar   }
1226e0b74bf9SHong Zhang 
12270298fd71SBarry 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);
1228d341cd04SHong 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() */
12290298fd71SBarry 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);
1230d341cd04SHong 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);
1231d341cd04SHong 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);
1232d341cd04SHong 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);
1233d341cd04SHong 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);
1234d341cd04SHong 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);
123559ac8732SStefano Zampini   if (mumps->id.ICNTL(19) <= 0 || mumps->id.ICNTL(19) > 3) { /* reset any schur data (if any) */
123659ac8732SStefano Zampini     ierr = MatMumpsResetSchur_Private(mumps);CHKERRQ(ierr);
123759ac8732SStefano Zampini   }
12384e34a73bSHong 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 */
1239d341cd04SHong 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 */
12409a2535b5SHong Zhang 
1241d341cd04SHong 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);
12420298fd71SBarry 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);
12430298fd71SBarry 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);
12449a2535b5SHong Zhang   if (mumps->id.ICNTL(24)) {
12459a2535b5SHong Zhang     mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */
1246d7ebd59bSHong Zhang   }
1247d7ebd59bSHong Zhang 
1248d341cd04SHong 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);
1249d341cd04SHong 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);
12502cd7d884SHong 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);
12510298fd71SBarry 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);
1252d341cd04SHong 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);
12530298fd71SBarry 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);
1254d341cd04SHong 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);
12554e34a73bSHong 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 */
12560298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),NULL);CHKERRQ(ierr);
1257dcd589f8SShri Abhyankar 
12580298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),NULL);CHKERRQ(ierr);
12590298fd71SBarry 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);
12600298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),NULL);CHKERRQ(ierr);
12610298fd71SBarry 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);
12620298fd71SBarry 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);
1263e5bb22a1SHong Zhang 
12640298fd71SBarry Smith   ierr = PetscOptionsString("-mat_mumps_ooc_tmpdir", "out of core directory", "None", mumps->id.ooc_tmpdir, mumps->id.ooc_tmpdir, 256, NULL);
1265b34f08ffSHong Zhang 
126616d797efSHong Zhang   ierr = PetscOptionsIntArray("-mat_mumps_view_info","request INFO local to each processor","",info,&ninfo,NULL);CHKERRQ(ierr);
1267b34f08ffSHong Zhang   if (ninfo) {
1268b34f08ffSHong Zhang     if (ninfo > 40) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"number of INFO %d must <= 40\n",ninfo);
1269b34f08ffSHong Zhang     ierr = PetscMalloc1(ninfo,&mumps->info);CHKERRQ(ierr);
1270b34f08ffSHong Zhang     mumps->ninfo = ninfo;
1271b34f08ffSHong Zhang     for (i=0; i<ninfo; i++) {
1272b34f08ffSHong Zhang       if (info[i] < 0 || info[i]>40) {
1273b34f08ffSHong Zhang         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index of INFO %d must between 1 and 40\n",ninfo);
1274b34f08ffSHong Zhang       } else {
1275b34f08ffSHong Zhang         mumps->info[i] = info[i];
1276b34f08ffSHong Zhang       }
1277b34f08ffSHong Zhang     }
1278b34f08ffSHong Zhang   }
1279b34f08ffSHong Zhang 
1280dcd589f8SShri Abhyankar   PetscOptionsEnd();
1281dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
1282dcd589f8SShri Abhyankar }
1283dcd589f8SShri Abhyankar 
1284dcd589f8SShri Abhyankar #undef __FUNCT__
1285dcd589f8SShri Abhyankar #define __FUNCT__ "PetscInitializeMUMPS"
1286f697e70eSHong Zhang PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps)
1287dcd589f8SShri Abhyankar {
1288dcd589f8SShri Abhyankar   PetscErrorCode ierr;
1289dcd589f8SShri Abhyankar 
1290dcd589f8SShri Abhyankar   PetscFunctionBegin;
1291ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)A), &mumps->myid);
1292ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&mumps->size);CHKERRQ(ierr);
1293ce94432eSBarry Smith   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(mumps->comm_mumps));CHKERRQ(ierr);
12942205254eSKarl Rupp 
1295f697e70eSHong Zhang   mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps);
1296f697e70eSHong Zhang 
1297f697e70eSHong Zhang   mumps->id.job = JOB_INIT;
1298f697e70eSHong Zhang   mumps->id.par = 1;  /* host participates factorizaton and solve */
1299f697e70eSHong Zhang   mumps->id.sym = mumps->sym;
13002907cef9SHong Zhang   PetscMUMPS_c(&mumps->id);
1301f697e70eSHong Zhang 
13020298fd71SBarry Smith   mumps->scat_rhs     = NULL;
13030298fd71SBarry Smith   mumps->scat_sol     = NULL;
13049a2535b5SHong Zhang 
130570544d5fSHong Zhang   /* set PETSc-MUMPS default options - override MUMPS default */
13069a2535b5SHong Zhang   mumps->id.ICNTL(3) = 0;
13079a2535b5SHong Zhang   mumps->id.ICNTL(4) = 0;
13089a2535b5SHong Zhang   if (mumps->size == 1) {
13099a2535b5SHong Zhang     mumps->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
13109a2535b5SHong Zhang   } else {
13119a2535b5SHong Zhang     mumps->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
13124e34a73bSHong Zhang     mumps->id.ICNTL(20) = 0;   /* rhs is in dense format */
131370544d5fSHong Zhang     mumps->id.ICNTL(21) = 1;   /* distributed solution */
13149a2535b5SHong Zhang   }
13156444a565SStefano Zampini 
13166444a565SStefano Zampini   /* schur */
13176444a565SStefano Zampini   mumps->id.size_schur      = 0;
13186444a565SStefano Zampini   mumps->id.listvar_schur   = NULL;
13196444a565SStefano Zampini   mumps->id.schur           = NULL;
1320b5fa320bSStefano Zampini   mumps->sizeredrhs         = 0;
1321b5fa320bSStefano Zampini   mumps->schur_pivots       = NULL;
1322b5fa320bSStefano Zampini   mumps->schur_work         = NULL;
132359ac8732SStefano Zampini   mumps->schur_sol          = NULL;
132459ac8732SStefano Zampini   mumps->schur_sizesol      = 0;
132559ac8732SStefano Zampini   mumps->schur_factored     = PETSC_FALSE;
132659ac8732SStefano Zampini   mumps->schur_inverted     = PETSC_FALSE;
1327dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
1328dcd589f8SShri Abhyankar }
1329dcd589f8SShri Abhyankar 
1330a5e57a09SHong 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 */
1331397b6df1SKris Buschelman #undef __FUNCT__
1332f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
13330481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1334b24902e0SBarry Smith {
1335a5e57a09SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
1336dcd589f8SShri Abhyankar   PetscErrorCode ierr;
133767877ebaSShri Abhyankar   Vec            b;
133867877ebaSShri Abhyankar   IS             is_iden;
133967877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1340397b6df1SKris Buschelman 
1341397b6df1SKris Buschelman   PetscFunctionBegin;
1342a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1343dcd589f8SShri Abhyankar 
13449a2535b5SHong Zhang   /* Set MUMPS options from the options database */
13459a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1346dcd589f8SShri Abhyankar 
1347a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1348dcd589f8SShri Abhyankar 
134967877ebaSShri Abhyankar   /* analysis phase */
135067877ebaSShri Abhyankar   /*----------------*/
1351a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1352a5e57a09SHong Zhang   mumps->id.n   = M;
1353a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
135467877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1355a5e57a09SHong Zhang     if (!mumps->myid) {
1356a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1357a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
1358940cd9d6SSatish Balay         mumps->id.a = (MumpsScalar*)mumps->val;
135967877ebaSShri Abhyankar       }
1360a5e57a09SHong Zhang       if (mumps->id.ICNTL(7) == 1) { /* use user-provide matrix ordering - assuming r = c ordering */
13615248a706SHong Zhang         /*
13625248a706SHong Zhang         PetscBool      flag;
13635248a706SHong Zhang         ierr = ISEqual(r,c,&flag);CHKERRQ(ierr);
13645248a706SHong Zhang         if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"row_perm != col_perm");
13655248a706SHong Zhang         ierr = ISView(r,PETSC_VIEWER_STDOUT_SELF);
13665248a706SHong Zhang          */
1367a5e57a09SHong Zhang         if (!mumps->myid) {
1368e0b74bf9SHong Zhang           const PetscInt *idx;
1369e0b74bf9SHong Zhang           PetscInt       i,*perm_in;
13702205254eSKarl Rupp 
1371785e854fSJed Brown           ierr = PetscMalloc1(M,&perm_in);CHKERRQ(ierr);
1372e0b74bf9SHong Zhang           ierr = ISGetIndices(r,&idx);CHKERRQ(ierr);
13732205254eSKarl Rupp 
1374a5e57a09SHong Zhang           mumps->id.perm_in = perm_in;
1375e0b74bf9SHong Zhang           for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */
1376e0b74bf9SHong Zhang           ierr = ISRestoreIndices(r,&idx);CHKERRQ(ierr);
1377e0b74bf9SHong Zhang         }
1378e0b74bf9SHong Zhang       }
137967877ebaSShri Abhyankar     }
138067877ebaSShri Abhyankar     break;
138167877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1382a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1383a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1384a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
1385940cd9d6SSatish Balay       mumps->id.a_loc = (MumpsScalar*)mumps->val;
138667877ebaSShri Abhyankar     }
138767877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1388a5e57a09SHong Zhang     if (!mumps->myid) {
13892cd7d884SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->rmap->N,&mumps->b_seq);CHKERRQ(ierr);
13902cd7d884SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->N,0,1,&is_iden);CHKERRQ(ierr);
139167877ebaSShri Abhyankar     } else {
1392a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
139367877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
139467877ebaSShri Abhyankar     }
13952a7a6963SBarry Smith     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1396a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
13976bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
13986bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
139967877ebaSShri Abhyankar     break;
140067877ebaSShri Abhyankar   }
1401a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
1402a5e57a09SHong 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));
140367877ebaSShri Abhyankar 
1404719d5645SBarry Smith   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1405dcd589f8SShri Abhyankar   F->ops->solve           = MatSolve_MUMPS;
140651d5961aSHong Zhang   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
14074e34a73bSHong Zhang   F->ops->matsolve        = MatMatSolve_MUMPS;
1408b24902e0SBarry Smith   PetscFunctionReturn(0);
1409b24902e0SBarry Smith }
1410b24902e0SBarry Smith 
1411450b117fSShri Abhyankar /* Note the Petsc r and c permutations are ignored */
1412450b117fSShri Abhyankar #undef __FUNCT__
1413450b117fSShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_BAIJMUMPS"
1414450b117fSShri Abhyankar PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1415450b117fSShri Abhyankar {
1416a5e57a09SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
1417dcd589f8SShri Abhyankar   PetscErrorCode ierr;
141867877ebaSShri Abhyankar   Vec            b;
141967877ebaSShri Abhyankar   IS             is_iden;
142067877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1421450b117fSShri Abhyankar 
1422450b117fSShri Abhyankar   PetscFunctionBegin;
1423a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1424dcd589f8SShri Abhyankar 
14259a2535b5SHong Zhang   /* Set MUMPS options from the options database */
14269a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1427dcd589f8SShri Abhyankar 
1428a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
142967877ebaSShri Abhyankar 
143067877ebaSShri Abhyankar   /* analysis phase */
143167877ebaSShri Abhyankar   /*----------------*/
1432a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1433a5e57a09SHong Zhang   mumps->id.n   = M;
1434a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
143567877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1436a5e57a09SHong Zhang     if (!mumps->myid) {
1437a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1438a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
1439940cd9d6SSatish Balay         mumps->id.a = (MumpsScalar*)mumps->val;
144067877ebaSShri Abhyankar       }
144167877ebaSShri Abhyankar     }
144267877ebaSShri Abhyankar     break;
144367877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1444a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1445a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1446a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
1447940cd9d6SSatish Balay       mumps->id.a_loc = (MumpsScalar*)mumps->val;
144867877ebaSShri Abhyankar     }
144967877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1450a5e57a09SHong Zhang     if (!mumps->myid) {
1451a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
145267877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
145367877ebaSShri Abhyankar     } else {
1454a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
145567877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
145667877ebaSShri Abhyankar     }
14572a7a6963SBarry Smith     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1458a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
14596bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
14606bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
146167877ebaSShri Abhyankar     break;
146267877ebaSShri Abhyankar   }
1463a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
1464a5e57a09SHong 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));
146567877ebaSShri Abhyankar 
1466450b117fSShri Abhyankar   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1467dcd589f8SShri Abhyankar   F->ops->solve           = MatSolve_MUMPS;
146851d5961aSHong Zhang   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1469450b117fSShri Abhyankar   PetscFunctionReturn(0);
1470450b117fSShri Abhyankar }
1471b24902e0SBarry Smith 
1472141f4205SHong Zhang /* Note the Petsc r permutation and factor info are ignored */
1473397b6df1SKris Buschelman #undef __FUNCT__
147467877ebaSShri Abhyankar #define __FUNCT__ "MatCholeskyFactorSymbolic_MUMPS"
147567877ebaSShri Abhyankar PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
1476b24902e0SBarry Smith {
1477a5e57a09SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
1478dcd589f8SShri Abhyankar   PetscErrorCode ierr;
147967877ebaSShri Abhyankar   Vec            b;
148067877ebaSShri Abhyankar   IS             is_iden;
148167877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1482397b6df1SKris Buschelman 
1483397b6df1SKris Buschelman   PetscFunctionBegin;
1484a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1485dcd589f8SShri Abhyankar 
14869a2535b5SHong Zhang   /* Set MUMPS options from the options database */
14879a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1488dcd589f8SShri Abhyankar 
1489a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1490dcd589f8SShri Abhyankar 
149167877ebaSShri Abhyankar   /* analysis phase */
149267877ebaSShri Abhyankar   /*----------------*/
1493a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1494a5e57a09SHong Zhang   mumps->id.n   = M;
1495a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
149667877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1497a5e57a09SHong Zhang     if (!mumps->myid) {
1498a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1499a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
1500940cd9d6SSatish Balay         mumps->id.a = (MumpsScalar*)mumps->val;
150167877ebaSShri Abhyankar       }
150267877ebaSShri Abhyankar     }
150367877ebaSShri Abhyankar     break;
150467877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1505a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1506a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1507a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
1508940cd9d6SSatish Balay       mumps->id.a_loc = (MumpsScalar*)mumps->val;
150967877ebaSShri Abhyankar     }
151067877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1511a5e57a09SHong Zhang     if (!mumps->myid) {
1512a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
151367877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
151467877ebaSShri Abhyankar     } else {
1515a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
151667877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
151767877ebaSShri Abhyankar     }
15182a7a6963SBarry Smith     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1519a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
15206bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
15216bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
152267877ebaSShri Abhyankar     break;
152367877ebaSShri Abhyankar   }
1524a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
1525a5e57a09SHong 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));
152667877ebaSShri Abhyankar 
15272792810eSHong Zhang   F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
1528dcd589f8SShri Abhyankar   F->ops->solve                 = MatSolve_MUMPS;
152951d5961aSHong Zhang   F->ops->solvetranspose        = MatSolve_MUMPS;
15304e34a73bSHong Zhang   F->ops->matsolve              = MatMatSolve_MUMPS;
15314e34a73bSHong Zhang #if defined(PETSC_USE_COMPLEX)
15320298fd71SBarry Smith   F->ops->getinertia = NULL;
15334e34a73bSHong Zhang #else
15344e34a73bSHong Zhang   F->ops->getinertia = MatGetInertia_SBAIJMUMPS;
1535db4efbfdSBarry Smith #endif
1536b24902e0SBarry Smith   PetscFunctionReturn(0);
1537b24902e0SBarry Smith }
1538b24902e0SBarry Smith 
1539397b6df1SKris Buschelman #undef __FUNCT__
154064e6c443SBarry Smith #define __FUNCT__ "MatView_MUMPS"
154164e6c443SBarry Smith PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
154274ed9c26SBarry Smith {
1543f6c57405SHong Zhang   PetscErrorCode    ierr;
154464e6c443SBarry Smith   PetscBool         iascii;
154564e6c443SBarry Smith   PetscViewerFormat format;
1546a5e57a09SHong Zhang   Mat_MUMPS         *mumps=(Mat_MUMPS*)A->spptr;
1547f6c57405SHong Zhang 
1548f6c57405SHong Zhang   PetscFunctionBegin;
154964e6c443SBarry Smith   /* check if matrix is mumps type */
155064e6c443SBarry Smith   if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0);
155164e6c443SBarry Smith 
1552251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
155364e6c443SBarry Smith   if (iascii) {
155464e6c443SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
155564e6c443SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO) {
155664e6c443SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
1557a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                   %d \n",mumps->id.sym);CHKERRQ(ierr);
1558a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):            %d \n",mumps->id.par);CHKERRQ(ierr);
1559a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):         %d \n",mumps->id.ICNTL(1));CHKERRQ(ierr);
1560a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg): %d \n",mumps->id.ICNTL(2));CHKERRQ(ierr);
1561a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):   %d \n",mumps->id.ICNTL(3));CHKERRQ(ierr);
1562a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):        %d \n",mumps->id.ICNTL(4));CHKERRQ(ierr);
1563a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):         %d \n",mumps->id.ICNTL(5));CHKERRQ(ierr);
1564a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):        %d \n",mumps->id.ICNTL(6));CHKERRQ(ierr);
1565a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (sequentia matrix ordering):%d \n",mumps->id.ICNTL(7));CHKERRQ(ierr);
1566a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scalling strategy):        %d \n",mumps->id.ICNTL(8));CHKERRQ(ierr);
1567a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements):  %d \n",mumps->id.ICNTL(10));CHKERRQ(ierr);
1568a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):          %d \n",mumps->id.ICNTL(11));CHKERRQ(ierr);
1569a5e57a09SHong Zhang       if (mumps->id.ICNTL(11)>0) {
1570a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(4) (inf norm of input mat):        %g\n",mumps->id.RINFOG(4));CHKERRQ(ierr);
1571a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(5) (inf norm of solution):         %g\n",mumps->id.RINFOG(5));CHKERRQ(ierr);
1572a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(6) (inf norm of residual):         %g\n",mumps->id.RINFOG(6));CHKERRQ(ierr);
1573a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8));CHKERRQ(ierr);
1574a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(9) (error estimate):               %g \n",mumps->id.RINFOG(9));CHKERRQ(ierr);
1575a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11));CHKERRQ(ierr);
1576f6c57405SHong Zhang       }
1577a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",mumps->id.ICNTL(12));CHKERRQ(ierr);
1578a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",mumps->id.ICNTL(13));CHKERRQ(ierr);
1579a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",mumps->id.ICNTL(14));CHKERRQ(ierr);
1580f6c57405SHong Zhang       /* ICNTL(15-17) not used */
1581a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",mumps->id.ICNTL(18));CHKERRQ(ierr);
1582a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Shur complement info):                       %d \n",mumps->id.ICNTL(19));CHKERRQ(ierr);
1583a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",mumps->id.ICNTL(20));CHKERRQ(ierr);
1584ca3dc52bSPierre Jolivet       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d \n",mumps->id.ICNTL(21));CHKERRQ(ierr);
1585a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",mumps->id.ICNTL(22));CHKERRQ(ierr);
1586a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",mumps->id.ICNTL(23));CHKERRQ(ierr);
1587c0165424SHong Zhang 
1588a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",mumps->id.ICNTL(24));CHKERRQ(ierr);
1589a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",mumps->id.ICNTL(25));CHKERRQ(ierr);
1590a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",mumps->id.ICNTL(26));CHKERRQ(ierr);
1591a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",mumps->id.ICNTL(27));CHKERRQ(ierr);
1592a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(28) (use parallel or sequential ordering):        %d \n",mumps->id.ICNTL(28));CHKERRQ(ierr);
1593a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(29) (parallel ordering):                          %d \n",mumps->id.ICNTL(29));CHKERRQ(ierr);
159442179a6aSHong Zhang 
1595a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(30) (user-specified set of entries in inv(A)):    %d \n",mumps->id.ICNTL(30));CHKERRQ(ierr);
1596a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(31) (factors is discarded in the solve phase):    %d \n",mumps->id.ICNTL(31));CHKERRQ(ierr);
1597a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(33) (compute determinant):                        %d \n",mumps->id.ICNTL(33));CHKERRQ(ierr);
1598f6c57405SHong Zhang 
1599a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",mumps->id.CNTL(1));CHKERRQ(ierr);
1600a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2));CHKERRQ(ierr);
1601ca3dc52bSPierre Jolivet       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",mumps->id.CNTL(3));CHKERRQ(ierr);
1602ca3dc52bSPierre Jolivet       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",mumps->id.CNTL(4));CHKERRQ(ierr);
1603a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",mumps->id.CNTL(5));CHKERRQ(ierr);
1604f6c57405SHong Zhang 
1605f6c57405SHong Zhang       /* infomation local to each processor */
160634ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);
16071575c14dSBarry Smith       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
1608a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %g \n",mumps->myid,mumps->id.RINFO(1));CHKERRQ(ierr);
160934ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
161034ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);
1611a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(2));CHKERRQ(ierr);
161234ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
161334ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);
1614a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(3));CHKERRQ(ierr);
161534ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1616f6c57405SHong Zhang 
161734ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr);
1618a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] %d \n",mumps->myid,mumps->id.INFO(15));CHKERRQ(ierr);
161934ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1620f6c57405SHong Zhang 
162134ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr);
1622a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(16));CHKERRQ(ierr);
162334ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1624f6c57405SHong Zhang 
162534ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);
1626a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(23));CHKERRQ(ierr);
162734ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1628b34f08ffSHong Zhang 
1629b34f08ffSHong Zhang       if (mumps->ninfo && mumps->ninfo <= 40){
1630b34f08ffSHong Zhang         PetscInt i;
1631b34f08ffSHong Zhang         for (i=0; i<mumps->ninfo; i++){
1632b34f08ffSHong Zhang           ierr = PetscViewerASCIIPrintf(viewer, "  INFO(%d): \n",mumps->info[i]);CHKERRQ(ierr);
1633b34f08ffSHong Zhang           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(mumps->info[i]));CHKERRQ(ierr);
1634b34f08ffSHong Zhang           ierr = PetscViewerFlush(viewer);
1635b34f08ffSHong Zhang         }
1636b34f08ffSHong Zhang       }
1637b34f08ffSHong Zhang 
1638b34f08ffSHong Zhang 
16391575c14dSBarry Smith       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
1640f6c57405SHong Zhang 
1641a5e57a09SHong Zhang       if (!mumps->myid) { /* information from the host */
1642a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",mumps->id.RINFOG(1));CHKERRQ(ierr);
1643a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",mumps->id.RINFOG(2));CHKERRQ(ierr);
1644a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",mumps->id.RINFOG(3));CHKERRQ(ierr);
1645a5e57a09SHong 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);
1646f6c57405SHong Zhang 
1647a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(3));CHKERRQ(ierr);
1648a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(4));CHKERRQ(ierr);
1649a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",mumps->id.INFOG(5));CHKERRQ(ierr);
1650a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",mumps->id.INFOG(6));CHKERRQ(ierr);
1651a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively use after analysis): %d \n",mumps->id.INFOG(7));CHKERRQ(ierr);
1652a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",mumps->id.INFOG(8));CHKERRQ(ierr);
1653a5e57a09SHong 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);
1654a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",mumps->id.INFOG(10));CHKERRQ(ierr);
1655a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",mumps->id.INFOG(11));CHKERRQ(ierr);
1656a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",mumps->id.INFOG(12));CHKERRQ(ierr);
1657a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",mumps->id.INFOG(13));CHKERRQ(ierr);
1658a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",mumps->id.INFOG(14));CHKERRQ(ierr);
1659a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",mumps->id.INFOG(15));CHKERRQ(ierr);
1660a5e57a09SHong 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);
1661a5e57a09SHong 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);
1662a5e57a09SHong 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);
1663a5e57a09SHong 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);
1664a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",mumps->id.INFOG(20));CHKERRQ(ierr);
1665a5e57a09SHong 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);
1666a5e57a09SHong 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);
1667a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",mumps->id.INFOG(23));CHKERRQ(ierr);
1668a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",mumps->id.INFOG(24));CHKERRQ(ierr);
1669a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",mumps->id.INFOG(25));CHKERRQ(ierr);
167040d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(28) (after factorization: number of null pivots encountered): %d\n",mumps->id.INFOG(28));CHKERRQ(ierr);
167140d435e3SHong 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);
167240d435e3SHong 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);
167340d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(32) (after analysis: type of analysis done): %d\n",mumps->id.INFOG(32));CHKERRQ(ierr);
167440d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(33) (value used for ICNTL(8)): %d\n",mumps->id.INFOG(33));CHKERRQ(ierr);
167540d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(34) (exponent of the determinant if determinant is requested): %d\n",mumps->id.INFOG(34));CHKERRQ(ierr);
1676f6c57405SHong Zhang       }
1677f6c57405SHong Zhang     }
1678cb828f0fSHong Zhang   }
1679f6c57405SHong Zhang   PetscFunctionReturn(0);
1680f6c57405SHong Zhang }
1681f6c57405SHong Zhang 
168235bd34faSBarry Smith #undef __FUNCT__
168335bd34faSBarry Smith #define __FUNCT__ "MatGetInfo_MUMPS"
168435bd34faSBarry Smith PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
168535bd34faSBarry Smith {
1686cb828f0fSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)A->spptr;
168735bd34faSBarry Smith 
168835bd34faSBarry Smith   PetscFunctionBegin;
168935bd34faSBarry Smith   info->block_size        = 1.0;
1690cb828f0fSHong Zhang   info->nz_allocated      = mumps->id.INFOG(20);
1691cb828f0fSHong Zhang   info->nz_used           = mumps->id.INFOG(20);
169235bd34faSBarry Smith   info->nz_unneeded       = 0.0;
169335bd34faSBarry Smith   info->assemblies        = 0.0;
169435bd34faSBarry Smith   info->mallocs           = 0.0;
169535bd34faSBarry Smith   info->memory            = 0.0;
169635bd34faSBarry Smith   info->fill_ratio_given  = 0;
169735bd34faSBarry Smith   info->fill_ratio_needed = 0;
169835bd34faSBarry Smith   info->factor_mallocs    = 0;
169935bd34faSBarry Smith   PetscFunctionReturn(0);
170035bd34faSBarry Smith }
170135bd34faSBarry Smith 
17025ccb76cbSHong Zhang /* -------------------------------------------------------------------------------------------*/
17035ccb76cbSHong Zhang #undef __FUNCT__
17048e7ba810SStefano Zampini #define __FUNCT__ "MatFactorSetSchurIS_MUMPS"
17058e7ba810SStefano Zampini PetscErrorCode MatFactorSetSchurIS_MUMPS(Mat F, IS is)
17066444a565SStefano Zampini {
17076444a565SStefano Zampini   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
17088e7ba810SStefano Zampini   const PetscInt *idxs;
17098e7ba810SStefano Zampini   PetscInt       size,i;
17106444a565SStefano Zampini   PetscErrorCode ierr;
17116444a565SStefano Zampini 
17126444a565SStefano Zampini   PetscFunctionBegin;
171359ac8732SStefano Zampini   if (mumps->size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MUMPS parallel Schur complements not yet supported from PETSc\n");
17148e7ba810SStefano Zampini   ierr = ISGetLocalSize(is,&size);CHKERRQ(ierr);
17156444a565SStefano Zampini   if (mumps->id.size_schur != size) {
17166444a565SStefano Zampini     ierr = PetscFree2(mumps->id.listvar_schur,mumps->id.schur);CHKERRQ(ierr);
17176444a565SStefano Zampini     mumps->id.size_schur = size;
17186444a565SStefano Zampini     mumps->id.schur_lld = size;
17196444a565SStefano Zampini     ierr = PetscMalloc2(size,&mumps->id.listvar_schur,size*size,&mumps->id.schur);CHKERRQ(ierr);
17206444a565SStefano Zampini   }
17218e7ba810SStefano Zampini   ierr = ISGetIndices(is,&idxs);CHKERRQ(ierr);
17226444a565SStefano Zampini   ierr = PetscMemcpy(mumps->id.listvar_schur,idxs,size*sizeof(PetscInt));CHKERRQ(ierr);
17238e7ba810SStefano Zampini   /* MUMPS expects Fortran style indices */
17248e7ba810SStefano Zampini   for (i=0;i<size;i++) mumps->id.listvar_schur[i]++;
17258e7ba810SStefano Zampini   ierr = ISRestoreIndices(is,&idxs);CHKERRQ(ierr);
1726*78a7176cSStefano Zampini   if (size) { /* turn on Schur switch if we the set of indices is not empty */
17276444a565SStefano Zampini     if (F->factortype == MAT_FACTOR_LU) {
172859ac8732SStefano Zampini       mumps->id.ICNTL(19) = 3; /* MUMPS returns full matrix */
17296444a565SStefano Zampini     } else {
173059ac8732SStefano Zampini       mumps->id.ICNTL(19) = 2; /* MUMPS returns lower triangular part */
17316444a565SStefano Zampini     }
173259ac8732SStefano Zampini     /* set a special value of ICNTL (not handled my MUMPS) to be used in the solve phase by PETSc */
1733b5fa320bSStefano Zampini     mumps->id.ICNTL(26) = -1;
1734*78a7176cSStefano Zampini   }
17356444a565SStefano Zampini   PetscFunctionReturn(0);
17366444a565SStefano Zampini }
17376444a565SStefano Zampini 
17386444a565SStefano Zampini /* -------------------------------------------------------------------------------------------*/
17396444a565SStefano Zampini #undef __FUNCT__
17405a05ddb0SStefano Zampini #define __FUNCT__ "MatFactorCreateSchurComplement_MUMPS"
17415a05ddb0SStefano Zampini PetscErrorCode MatFactorCreateSchurComplement_MUMPS(Mat F,Mat* S)
17426444a565SStefano Zampini {
17436444a565SStefano Zampini   Mat            St;
17446444a565SStefano Zampini   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
17456444a565SStefano Zampini   PetscScalar    *array;
17466444a565SStefano Zampini #if defined(PETSC_USE_COMPLEX)
17478ac429a0SStefano Zampini   PetscScalar    im = PetscSqrtScalar((PetscScalar)-1.0);
17486444a565SStefano Zampini #endif
17496444a565SStefano Zampini   PetscErrorCode ierr;
17506444a565SStefano Zampini 
17516444a565SStefano Zampini   PetscFunctionBegin;
17525a05ddb0SStefano 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");
17535a05ddb0SStefano Zampini   else if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before");
17546f3cc6f9SBarry Smith 
17556444a565SStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)F),&St);CHKERRQ(ierr);
17566444a565SStefano Zampini   ierr = MatSetSizes(St,PETSC_DECIDE,PETSC_DECIDE,mumps->id.size_schur,mumps->id.size_schur);CHKERRQ(ierr);
17576444a565SStefano Zampini   ierr = MatSetType(St,MATDENSE);CHKERRQ(ierr);
17586444a565SStefano Zampini   ierr = MatSetUp(St);CHKERRQ(ierr);
17596444a565SStefano Zampini   ierr = MatDenseGetArray(St,&array);CHKERRQ(ierr);
176059ac8732SStefano Zampini   if (!mumps->sym) { /* MUMPS always return a full matrix */
17616444a565SStefano Zampini     if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
17626444a565SStefano Zampini       PetscInt i,j,N=mumps->id.size_schur;
17636444a565SStefano Zampini       for (i=0;i<N;i++) {
17646444a565SStefano Zampini         for (j=0;j<N;j++) {
17656444a565SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
17666444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j];
17676444a565SStefano Zampini #else
17686444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
17696444a565SStefano Zampini #endif
17706444a565SStefano Zampini           array[j*N+i] = val;
17716444a565SStefano Zampini         }
17726444a565SStefano Zampini       }
17736444a565SStefano Zampini     } else { /* stored by columns */
17746444a565SStefano Zampini       ierr = PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
17756444a565SStefano Zampini     }
17766444a565SStefano Zampini   } else { /* either full or lower-triangular (not packed) */
17776444a565SStefano Zampini     if (mumps->id.ICNTL(19) == 2) { /* lower triangular stored by columns */
17786444a565SStefano Zampini       PetscInt i,j,N=mumps->id.size_schur;
17796444a565SStefano Zampini       for (i=0;i<N;i++) {
17806444a565SStefano Zampini         for (j=i;j<N;j++) {
17816444a565SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
17826444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j];
17836444a565SStefano Zampini #else
17846444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
17856444a565SStefano Zampini #endif
17866444a565SStefano Zampini           array[i*N+j] = val;
17876444a565SStefano Zampini           array[j*N+i] = val;
17886444a565SStefano Zampini         }
17896444a565SStefano Zampini       }
17906444a565SStefano Zampini     } else if (mumps->id.ICNTL(19) == 3) { /* full matrix */
17916444a565SStefano Zampini       ierr = PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
17926444a565SStefano Zampini     } else { /* ICNTL(19) == 1 lower triangular stored by rows */
17936444a565SStefano Zampini       PetscInt i,j,N=mumps->id.size_schur;
17946444a565SStefano Zampini       for (i=0;i<N;i++) {
17956444a565SStefano Zampini         for (j=0;j<i+1;j++) {
17966444a565SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
17976444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j];
17986444a565SStefano Zampini #else
17996444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
18006444a565SStefano Zampini #endif
18016444a565SStefano Zampini           array[i*N+j] = val;
18026444a565SStefano Zampini           array[j*N+i] = val;
18036444a565SStefano Zampini         }
18046444a565SStefano Zampini       }
18056444a565SStefano Zampini     }
18066444a565SStefano Zampini   }
18076444a565SStefano Zampini   ierr = MatDenseRestoreArray(St,&array);CHKERRQ(ierr);
18086444a565SStefano Zampini   *S = St;
18096444a565SStefano Zampini   PetscFunctionReturn(0);
18106444a565SStefano Zampini }
18116444a565SStefano Zampini 
181259ac8732SStefano Zampini /* -------------------------------------------------------------------------------------------*/
181359ac8732SStefano Zampini #undef __FUNCT__
18145a05ddb0SStefano Zampini #define __FUNCT__ "MatFactorGetSchurComplement_MUMPS"
18155a05ddb0SStefano Zampini PetscErrorCode MatFactorGetSchurComplement_MUMPS(Mat F,Mat* S)
181659ac8732SStefano Zampini {
181759ac8732SStefano Zampini   Mat            St;
181859ac8732SStefano Zampini   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
181959ac8732SStefano Zampini   PetscErrorCode ierr;
182059ac8732SStefano Zampini 
182159ac8732SStefano Zampini   PetscFunctionBegin;
18225a05ddb0SStefano 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");
18235a05ddb0SStefano Zampini   else if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before");
18246f3cc6f9SBarry Smith 
18257bc32933SStefano 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 */
182659ac8732SStefano Zampini   ierr = MatCreateSeqDense(PetscObjectComm((PetscObject)F),mumps->id.size_schur,mumps->id.size_schur,(PetscScalar*)mumps->id.schur,&St);CHKERRQ(ierr);
182759ac8732SStefano Zampini   *S = St;
182859ac8732SStefano Zampini   PetscFunctionReturn(0);
182959ac8732SStefano Zampini }
183059ac8732SStefano Zampini 
183159ac8732SStefano Zampini /* -------------------------------------------------------------------------------------------*/
183259ac8732SStefano Zampini #undef __FUNCT__
18335a05ddb0SStefano Zampini #define __FUNCT__ "MatFactorInvertSchurComplement_MUMPS"
18345a05ddb0SStefano Zampini PetscErrorCode MatFactorInvertSchurComplement_MUMPS(Mat F)
183559ac8732SStefano Zampini {
183659ac8732SStefano Zampini   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
183759ac8732SStefano Zampini   PetscErrorCode ierr;
183859ac8732SStefano Zampini 
183959ac8732SStefano Zampini   PetscFunctionBegin;
184059ac8732SStefano Zampini   if (!mumps->id.ICNTL(19)) { /* do nothing */
184159ac8732SStefano Zampini     PetscFunctionReturn(0);
184259ac8732SStefano Zampini   }
18435a05ddb0SStefano Zampini   if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before");
184459ac8732SStefano Zampini   ierr = MatMumpsInvertSchur_Private(mumps);CHKERRQ(ierr);
184559ac8732SStefano Zampini   PetscFunctionReturn(0);
184659ac8732SStefano Zampini }
184759ac8732SStefano Zampini 
18486444a565SStefano Zampini /* -------------------------------------------------------------------------------------------*/
18496444a565SStefano Zampini #undef __FUNCT__
18505a05ddb0SStefano Zampini #define __FUNCT__ "MatFactorSolveSchurComplement_MUMPS"
18515a05ddb0SStefano Zampini PetscErrorCode MatFactorSolveSchurComplement_MUMPS(Mat F, Vec rhs, Vec sol)
1852e807eca7SStefano Zampini {
1853e807eca7SStefano Zampini   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
1854e807eca7SStefano Zampini   MumpsScalar    *orhs;
1855e807eca7SStefano Zampini   PetscScalar    *osol,*nrhs,*nsol;
18569a3a5937SStefano Zampini   PetscInt       orhs_size,osol_size,olrhs_size;
1857e807eca7SStefano Zampini   PetscErrorCode ierr;
1858e807eca7SStefano Zampini 
1859e807eca7SStefano Zampini   PetscFunctionBegin;
18605a05ddb0SStefano 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");
18615a05ddb0SStefano Zampini   if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before");
18626f3cc6f9SBarry Smith 
1863e807eca7SStefano Zampini   /* swap pointers */
1864e807eca7SStefano Zampini   orhs = mumps->id.redrhs;
18659a3a5937SStefano Zampini   olrhs_size = mumps->id.lredrhs;
1866a12f35bfSStefano Zampini   orhs_size = mumps->sizeredrhs;
1867e807eca7SStefano Zampini   osol = mumps->schur_sol;
1868a12f35bfSStefano Zampini   osol_size = mumps->schur_sizesol;
1869e807eca7SStefano Zampini   ierr = VecGetArray(rhs,&nrhs);CHKERRQ(ierr);
1870e807eca7SStefano Zampini   ierr = VecGetArray(sol,&nsol);CHKERRQ(ierr);
1871e807eca7SStefano Zampini   mumps->id.redrhs = (MumpsScalar*)nrhs;
1872a12f35bfSStefano Zampini   ierr = VecGetLocalSize(rhs,&mumps->sizeredrhs);CHKERRQ(ierr);
18739a3a5937SStefano Zampini   mumps->id.lredrhs = mumps->sizeredrhs;
1874e807eca7SStefano Zampini   mumps->schur_sol = nsol;
1875a12f35bfSStefano Zampini   ierr = VecGetLocalSize(sol,&mumps->schur_sizesol);CHKERRQ(ierr);
1876a12f35bfSStefano Zampini 
1877e807eca7SStefano Zampini   /* solve Schur complement */
1878e807eca7SStefano Zampini   mumps->id.nrhs = 1;
1879e807eca7SStefano Zampini   ierr = MatMumpsSolveSchur_Private(mumps,PETSC_FALSE);CHKERRQ(ierr);
1880e807eca7SStefano Zampini   /* restore pointers */
1881e807eca7SStefano Zampini   ierr = VecRestoreArray(rhs,&nrhs);CHKERRQ(ierr);
1882e807eca7SStefano Zampini   ierr = VecRestoreArray(sol,&nsol);CHKERRQ(ierr);
1883e807eca7SStefano Zampini   mumps->id.redrhs = orhs;
18849a3a5937SStefano Zampini   mumps->id.lredrhs = olrhs_size;
1885a12f35bfSStefano Zampini   mumps->sizeredrhs = orhs_size;
1886e807eca7SStefano Zampini   mumps->schur_sol = osol;
1887a12f35bfSStefano Zampini   mumps->schur_sizesol = osol_size;
1888e807eca7SStefano Zampini   PetscFunctionReturn(0);
1889e807eca7SStefano Zampini }
1890e807eca7SStefano Zampini 
18917404bcfbSStefano Zampini /* -------------------------------------------------------------------------------------------*/
18927404bcfbSStefano Zampini #undef __FUNCT__
18935a05ddb0SStefano Zampini #define __FUNCT__ "MatFactorSolveSchurComplementTranspose_MUMPS"
18945a05ddb0SStefano Zampini PetscErrorCode MatFactorSolveSchurComplementTranspose_MUMPS(Mat F, Vec rhs, Vec sol)
18957404bcfbSStefano Zampini {
18967404bcfbSStefano Zampini   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
18977404bcfbSStefano Zampini   MumpsScalar    *orhs;
18987404bcfbSStefano Zampini   PetscScalar    *osol,*nrhs,*nsol;
1899a12f35bfSStefano Zampini   PetscInt       orhs_size,osol_size;
19007404bcfbSStefano Zampini   PetscErrorCode ierr;
19017404bcfbSStefano Zampini 
19027404bcfbSStefano Zampini   PetscFunctionBegin;
19035a05ddb0SStefano 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");
19045a05ddb0SStefano Zampini   else if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before");
19056f3cc6f9SBarry Smith 
19067404bcfbSStefano Zampini   /* swap pointers */
19077404bcfbSStefano Zampini   orhs = mumps->id.redrhs;
1908a12f35bfSStefano Zampini   orhs_size = mumps->sizeredrhs;
19097404bcfbSStefano Zampini   osol = mumps->schur_sol;
1910a12f35bfSStefano Zampini   osol_size = mumps->schur_sizesol;
19117404bcfbSStefano Zampini   ierr = VecGetArray(rhs,&nrhs);CHKERRQ(ierr);
19127404bcfbSStefano Zampini   ierr = VecGetArray(sol,&nsol);CHKERRQ(ierr);
19137404bcfbSStefano Zampini   mumps->id.redrhs = (MumpsScalar*)nrhs;
1914a12f35bfSStefano Zampini   ierr = VecGetLocalSize(rhs,&mumps->sizeredrhs);CHKERRQ(ierr);
19157404bcfbSStefano Zampini   mumps->schur_sol = nsol;
1916a12f35bfSStefano Zampini   ierr = VecGetLocalSize(sol,&mumps->schur_sizesol);CHKERRQ(ierr);
1917a12f35bfSStefano Zampini 
19187404bcfbSStefano Zampini   /* solve Schur complement */
19197404bcfbSStefano Zampini   mumps->id.nrhs = 1;
19207404bcfbSStefano Zampini   mumps->id.ICNTL(9) = 0;
19217404bcfbSStefano Zampini   ierr = MatMumpsSolveSchur_Private(mumps,PETSC_FALSE);CHKERRQ(ierr);
19227404bcfbSStefano Zampini   mumps->id.ICNTL(9) = 1;
19237404bcfbSStefano Zampini   /* restore pointers */
19247404bcfbSStefano Zampini   ierr = VecRestoreArray(rhs,&nrhs);CHKERRQ(ierr);
19257404bcfbSStefano Zampini   ierr = VecRestoreArray(sol,&nsol);CHKERRQ(ierr);
19267404bcfbSStefano Zampini   mumps->id.redrhs = orhs;
1927a12f35bfSStefano Zampini   mumps->sizeredrhs = orhs_size;
19287404bcfbSStefano Zampini   mumps->schur_sol = osol;
1929a12f35bfSStefano Zampini   mumps->schur_sizesol = osol_size;
19307404bcfbSStefano Zampini   PetscFunctionReturn(0);
19317404bcfbSStefano Zampini }
19327404bcfbSStefano Zampini 
1933e807eca7SStefano Zampini /* -------------------------------------------------------------------------------------------*/
1934e807eca7SStefano Zampini #undef __FUNCT__
19355ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl_MUMPS"
19365ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival)
19375ccb76cbSHong Zhang {
1938a5e57a09SHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
19395ccb76cbSHong Zhang 
19405ccb76cbSHong Zhang   PetscFunctionBegin;
1941a5e57a09SHong Zhang   mumps->id.ICNTL(icntl) = ival;
19425ccb76cbSHong Zhang   PetscFunctionReturn(0);
19435ccb76cbSHong Zhang }
19445ccb76cbSHong Zhang 
19455ccb76cbSHong Zhang #undef __FUNCT__
1946bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetIcntl_MUMPS"
1947bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt *ival)
1948bc6112feSHong Zhang {
1949bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1950bc6112feSHong Zhang 
1951bc6112feSHong Zhang   PetscFunctionBegin;
1952bc6112feSHong Zhang   *ival = mumps->id.ICNTL(icntl);
1953bc6112feSHong Zhang   PetscFunctionReturn(0);
1954bc6112feSHong Zhang }
1955bc6112feSHong Zhang 
1956bc6112feSHong Zhang #undef __FUNCT__
19575ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl"
19585ccb76cbSHong Zhang /*@
19595ccb76cbSHong Zhang   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()
19605ccb76cbSHong Zhang 
19615ccb76cbSHong Zhang    Logically Collective on Mat
19625ccb76cbSHong Zhang 
19635ccb76cbSHong Zhang    Input Parameters:
19645ccb76cbSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
19655ccb76cbSHong Zhang .  icntl - index of MUMPS parameter array ICNTL()
19665ccb76cbSHong Zhang -  ival - value of MUMPS ICNTL(icntl)
19675ccb76cbSHong Zhang 
19685ccb76cbSHong Zhang   Options Database:
19695ccb76cbSHong Zhang .   -mat_mumps_icntl_<icntl> <ival>
19705ccb76cbSHong Zhang 
19715ccb76cbSHong Zhang    Level: beginner
19725ccb76cbSHong Zhang 
19735ccb76cbSHong Zhang    References: MUMPS Users' Guide
19745ccb76cbSHong Zhang 
19755ccb76cbSHong Zhang .seealso: MatGetFactor()
19765ccb76cbSHong Zhang @*/
19775ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival)
19785ccb76cbSHong Zhang {
19795ccb76cbSHong Zhang   PetscErrorCode ierr;
19805ccb76cbSHong Zhang 
19815ccb76cbSHong Zhang   PetscFunctionBegin;
19825ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
19835ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,ival,3);
19845ccb76cbSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr);
19855ccb76cbSHong Zhang   PetscFunctionReturn(0);
19865ccb76cbSHong Zhang }
19875ccb76cbSHong Zhang 
1988bc6112feSHong Zhang #undef __FUNCT__
1989bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetIcntl"
1990a21f80fcSHong Zhang /*@
1991a21f80fcSHong Zhang   MatMumpsGetIcntl - Get MUMPS parameter ICNTL()
1992a21f80fcSHong Zhang 
1993a21f80fcSHong Zhang    Logically Collective on Mat
1994a21f80fcSHong Zhang 
1995a21f80fcSHong Zhang    Input Parameters:
1996a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1997a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array ICNTL()
1998a21f80fcSHong Zhang 
1999a21f80fcSHong Zhang   Output Parameter:
2000a21f80fcSHong Zhang .  ival - value of MUMPS ICNTL(icntl)
2001a21f80fcSHong Zhang 
2002a21f80fcSHong Zhang    Level: beginner
2003a21f80fcSHong Zhang 
2004a21f80fcSHong Zhang    References: MUMPS Users' Guide
2005a21f80fcSHong Zhang 
2006a21f80fcSHong Zhang .seealso: MatGetFactor()
2007a21f80fcSHong Zhang @*/
2008bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl(Mat F,PetscInt icntl,PetscInt *ival)
2009bc6112feSHong Zhang {
2010bc6112feSHong Zhang   PetscErrorCode ierr;
2011bc6112feSHong Zhang 
2012bc6112feSHong Zhang   PetscFunctionBegin;
2013bc6112feSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
2014bc6112feSHong Zhang   PetscValidIntPointer(ival,3);
2015bc6112feSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetIcntl_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
2016bc6112feSHong Zhang   PetscFunctionReturn(0);
2017bc6112feSHong Zhang }
2018bc6112feSHong Zhang 
20198928b65cSHong Zhang /* -------------------------------------------------------------------------------------------*/
20208928b65cSHong Zhang #undef __FUNCT__
20218928b65cSHong Zhang #define __FUNCT__ "MatMumpsSetCntl_MUMPS"
20228928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal val)
20238928b65cSHong Zhang {
20248928b65cSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
20258928b65cSHong Zhang 
20268928b65cSHong Zhang   PetscFunctionBegin;
20278928b65cSHong Zhang   mumps->id.CNTL(icntl) = val;
20288928b65cSHong Zhang   PetscFunctionReturn(0);
20298928b65cSHong Zhang }
20308928b65cSHong Zhang 
20318928b65cSHong Zhang #undef __FUNCT__
2032bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetCntl_MUMPS"
2033bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal *val)
2034bc6112feSHong Zhang {
2035bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
2036bc6112feSHong Zhang 
2037bc6112feSHong Zhang   PetscFunctionBegin;
2038bc6112feSHong Zhang   *val = mumps->id.CNTL(icntl);
2039bc6112feSHong Zhang   PetscFunctionReturn(0);
2040bc6112feSHong Zhang }
2041bc6112feSHong Zhang 
2042bc6112feSHong Zhang #undef __FUNCT__
20438928b65cSHong Zhang #define __FUNCT__ "MatMumpsSetCntl"
20448928b65cSHong Zhang /*@
20458928b65cSHong Zhang   MatMumpsSetCntl - Set MUMPS parameter CNTL()
20468928b65cSHong Zhang 
20478928b65cSHong Zhang    Logically Collective on Mat
20488928b65cSHong Zhang 
20498928b65cSHong Zhang    Input Parameters:
20508928b65cSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
20518928b65cSHong Zhang .  icntl - index of MUMPS parameter array CNTL()
20528928b65cSHong Zhang -  val - value of MUMPS CNTL(icntl)
20538928b65cSHong Zhang 
20548928b65cSHong Zhang   Options Database:
20558928b65cSHong Zhang .   -mat_mumps_cntl_<icntl> <val>
20568928b65cSHong Zhang 
20578928b65cSHong Zhang    Level: beginner
20588928b65cSHong Zhang 
20598928b65cSHong Zhang    References: MUMPS Users' Guide
20608928b65cSHong Zhang 
20618928b65cSHong Zhang .seealso: MatGetFactor()
20628928b65cSHong Zhang @*/
20638928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl(Mat F,PetscInt icntl,PetscReal val)
20648928b65cSHong Zhang {
20658928b65cSHong Zhang   PetscErrorCode ierr;
20668928b65cSHong Zhang 
20678928b65cSHong Zhang   PetscFunctionBegin;
20688928b65cSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
2069bc6112feSHong Zhang   PetscValidLogicalCollectiveReal(F,val,3);
20708928b65cSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsSetCntl_C",(Mat,PetscInt,PetscReal),(F,icntl,val));CHKERRQ(ierr);
20718928b65cSHong Zhang   PetscFunctionReturn(0);
20728928b65cSHong Zhang }
20738928b65cSHong Zhang 
2074bc6112feSHong Zhang #undef __FUNCT__
2075bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetCntl"
2076a21f80fcSHong Zhang /*@
2077a21f80fcSHong Zhang   MatMumpsGetCntl - Get MUMPS parameter CNTL()
2078a21f80fcSHong Zhang 
2079a21f80fcSHong Zhang    Logically Collective on Mat
2080a21f80fcSHong Zhang 
2081a21f80fcSHong Zhang    Input Parameters:
2082a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2083a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array CNTL()
2084a21f80fcSHong Zhang 
2085a21f80fcSHong Zhang   Output Parameter:
2086a21f80fcSHong Zhang .  val - value of MUMPS CNTL(icntl)
2087a21f80fcSHong Zhang 
2088a21f80fcSHong Zhang    Level: beginner
2089a21f80fcSHong Zhang 
2090a21f80fcSHong Zhang    References: MUMPS Users' Guide
2091a21f80fcSHong Zhang 
2092a21f80fcSHong Zhang .seealso: MatGetFactor()
2093a21f80fcSHong Zhang @*/
2094bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl(Mat F,PetscInt icntl,PetscReal *val)
2095bc6112feSHong Zhang {
2096bc6112feSHong Zhang   PetscErrorCode ierr;
2097bc6112feSHong Zhang 
2098bc6112feSHong Zhang   PetscFunctionBegin;
2099bc6112feSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
2100bc6112feSHong Zhang   PetscValidRealPointer(val,3);
2101bc6112feSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetCntl_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
2102bc6112feSHong Zhang   PetscFunctionReturn(0);
2103bc6112feSHong Zhang }
2104bc6112feSHong Zhang 
2105bc6112feSHong Zhang #undef __FUNCT__
2106ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfo_MUMPS"
2107ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F,PetscInt icntl,PetscInt *info)
2108bc6112feSHong Zhang {
2109bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
2110bc6112feSHong Zhang 
2111bc6112feSHong Zhang   PetscFunctionBegin;
2112bc6112feSHong Zhang   *info = mumps->id.INFO(icntl);
2113bc6112feSHong Zhang   PetscFunctionReturn(0);
2114bc6112feSHong Zhang }
2115bc6112feSHong Zhang 
2116bc6112feSHong Zhang #undef __FUNCT__
2117ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfog_MUMPS"
2118ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F,PetscInt icntl,PetscInt *infog)
2119bc6112feSHong Zhang {
2120bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
2121bc6112feSHong Zhang 
2122bc6112feSHong Zhang   PetscFunctionBegin;
2123bc6112feSHong Zhang   *infog = mumps->id.INFOG(icntl);
2124bc6112feSHong Zhang   PetscFunctionReturn(0);
2125bc6112feSHong Zhang }
2126bc6112feSHong Zhang 
2127bc6112feSHong Zhang #undef __FUNCT__
2128ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfo_MUMPS"
2129ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfo)
2130bc6112feSHong Zhang {
2131bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
2132bc6112feSHong Zhang 
2133bc6112feSHong Zhang   PetscFunctionBegin;
2134bc6112feSHong Zhang   *rinfo = mumps->id.RINFO(icntl);
2135bc6112feSHong Zhang   PetscFunctionReturn(0);
2136bc6112feSHong Zhang }
2137bc6112feSHong Zhang 
2138bc6112feSHong Zhang #undef __FUNCT__
2139ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfog_MUMPS"
2140ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfog)
2141bc6112feSHong Zhang {
2142bc6112feSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
2143bc6112feSHong Zhang 
2144bc6112feSHong Zhang   PetscFunctionBegin;
2145bc6112feSHong Zhang   *rinfog = mumps->id.RINFOG(icntl);
2146bc6112feSHong Zhang   PetscFunctionReturn(0);
2147bc6112feSHong Zhang }
2148bc6112feSHong Zhang 
2149bc6112feSHong Zhang #undef __FUNCT__
2150ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfo"
2151a21f80fcSHong Zhang /*@
2152a21f80fcSHong Zhang   MatMumpsGetInfo - Get MUMPS parameter INFO()
2153a21f80fcSHong Zhang 
2154a21f80fcSHong Zhang    Logically Collective on Mat
2155a21f80fcSHong Zhang 
2156a21f80fcSHong Zhang    Input Parameters:
2157a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2158a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array INFO()
2159a21f80fcSHong Zhang 
2160a21f80fcSHong Zhang   Output Parameter:
2161a21f80fcSHong Zhang .  ival - value of MUMPS INFO(icntl)
2162a21f80fcSHong Zhang 
2163a21f80fcSHong Zhang    Level: beginner
2164a21f80fcSHong Zhang 
2165a21f80fcSHong Zhang    References: MUMPS Users' Guide
2166a21f80fcSHong Zhang 
2167a21f80fcSHong Zhang .seealso: MatGetFactor()
2168a21f80fcSHong Zhang @*/
2169ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo(Mat F,PetscInt icntl,PetscInt *ival)
2170bc6112feSHong Zhang {
2171bc6112feSHong Zhang   PetscErrorCode ierr;
2172bc6112feSHong Zhang 
2173bc6112feSHong Zhang   PetscFunctionBegin;
2174ca810319SHong Zhang   PetscValidIntPointer(ival,3);
2175ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetInfo_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
2176bc6112feSHong Zhang   PetscFunctionReturn(0);
2177bc6112feSHong Zhang }
2178bc6112feSHong Zhang 
2179bc6112feSHong Zhang #undef __FUNCT__
2180ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfog"
2181a21f80fcSHong Zhang /*@
2182a21f80fcSHong Zhang   MatMumpsGetInfog - Get MUMPS parameter INFOG()
2183a21f80fcSHong Zhang 
2184a21f80fcSHong Zhang    Logically Collective on Mat
2185a21f80fcSHong Zhang 
2186a21f80fcSHong Zhang    Input Parameters:
2187a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2188a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array INFOG()
2189a21f80fcSHong Zhang 
2190a21f80fcSHong Zhang   Output Parameter:
2191a21f80fcSHong Zhang .  ival - value of MUMPS INFOG(icntl)
2192a21f80fcSHong Zhang 
2193a21f80fcSHong Zhang    Level: beginner
2194a21f80fcSHong Zhang 
2195a21f80fcSHong Zhang    References: MUMPS Users' Guide
2196a21f80fcSHong Zhang 
2197a21f80fcSHong Zhang .seealso: MatGetFactor()
2198a21f80fcSHong Zhang @*/
2199ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog(Mat F,PetscInt icntl,PetscInt *ival)
2200bc6112feSHong Zhang {
2201bc6112feSHong Zhang   PetscErrorCode ierr;
2202bc6112feSHong Zhang 
2203bc6112feSHong Zhang   PetscFunctionBegin;
2204ca810319SHong Zhang   PetscValidIntPointer(ival,3);
2205ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetInfog_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
2206bc6112feSHong Zhang   PetscFunctionReturn(0);
2207bc6112feSHong Zhang }
2208bc6112feSHong Zhang 
2209bc6112feSHong Zhang #undef __FUNCT__
2210ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfo"
2211a21f80fcSHong Zhang /*@
2212a21f80fcSHong Zhang   MatMumpsGetRinfo - Get MUMPS parameter RINFO()
2213a21f80fcSHong Zhang 
2214a21f80fcSHong Zhang    Logically Collective on Mat
2215a21f80fcSHong Zhang 
2216a21f80fcSHong Zhang    Input Parameters:
2217a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2218a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array RINFO()
2219a21f80fcSHong Zhang 
2220a21f80fcSHong Zhang   Output Parameter:
2221a21f80fcSHong Zhang .  val - value of MUMPS RINFO(icntl)
2222a21f80fcSHong Zhang 
2223a21f80fcSHong Zhang    Level: beginner
2224a21f80fcSHong Zhang 
2225a21f80fcSHong Zhang    References: MUMPS Users' Guide
2226a21f80fcSHong Zhang 
2227a21f80fcSHong Zhang .seealso: MatGetFactor()
2228a21f80fcSHong Zhang @*/
2229ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo(Mat F,PetscInt icntl,PetscReal *val)
2230bc6112feSHong Zhang {
2231bc6112feSHong Zhang   PetscErrorCode ierr;
2232bc6112feSHong Zhang 
2233bc6112feSHong Zhang   PetscFunctionBegin;
2234bc6112feSHong Zhang   PetscValidRealPointer(val,3);
2235ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetRinfo_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
2236bc6112feSHong Zhang   PetscFunctionReturn(0);
2237bc6112feSHong Zhang }
2238bc6112feSHong Zhang 
2239bc6112feSHong Zhang #undef __FUNCT__
2240ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfog"
2241a21f80fcSHong Zhang /*@
2242a21f80fcSHong Zhang   MatMumpsGetRinfog - Get MUMPS parameter RINFOG()
2243a21f80fcSHong Zhang 
2244a21f80fcSHong Zhang    Logically Collective on Mat
2245a21f80fcSHong Zhang 
2246a21f80fcSHong Zhang    Input Parameters:
2247a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2248a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array RINFOG()
2249a21f80fcSHong Zhang 
2250a21f80fcSHong Zhang   Output Parameter:
2251a21f80fcSHong Zhang .  val - value of MUMPS RINFOG(icntl)
2252a21f80fcSHong Zhang 
2253a21f80fcSHong Zhang    Level: beginner
2254a21f80fcSHong Zhang 
2255a21f80fcSHong Zhang    References: MUMPS Users' Guide
2256a21f80fcSHong Zhang 
2257a21f80fcSHong Zhang .seealso: MatGetFactor()
2258a21f80fcSHong Zhang @*/
2259ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog(Mat F,PetscInt icntl,PetscReal *val)
2260bc6112feSHong Zhang {
2261bc6112feSHong Zhang   PetscErrorCode ierr;
2262bc6112feSHong Zhang 
2263bc6112feSHong Zhang   PetscFunctionBegin;
2264bc6112feSHong Zhang   PetscValidRealPointer(val,3);
2265ca810319SHong Zhang   ierr = PetscTryMethod(F,"MatMumpsGetRinfog_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
2266bc6112feSHong Zhang   PetscFunctionReturn(0);
2267bc6112feSHong Zhang }
2268bc6112feSHong Zhang 
226924b6179bSKris Buschelman /*MC
22702692d6eeSBarry Smith   MATSOLVERMUMPS -  A matrix type providing direct solvers (LU and Cholesky) for
227124b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
227224b6179bSKris Buschelman 
227341c8de11SBarry Smith   Works with MATAIJ and MATSBAIJ matrices
227424b6179bSKris Buschelman 
2275c2b89b5dSBarry Smith   Use ./configure --download-mumps --download-scalapack --download-parmetis --download-metis --download-ptscotch  to have PETSc installed with MUMPS
2276c2b89b5dSBarry Smith 
2277c2b89b5dSBarry Smith   Use -pc_type cholesky or lu -pc_factor_mat_solver_package mumps to us this direct solver
2278c2b89b5dSBarry Smith 
227924b6179bSKris Buschelman   Options Database Keys:
22804e34a73bSHong Zhang +  -mat_mumps_icntl_1 <6>: ICNTL(1): output stream for error messages (None)
22814e34a73bSHong Zhang .  -mat_mumps_icntl_2 <0>: ICNTL(2): output stream for diagnostic printing, statistics, and warning (None)
22824e34a73bSHong Zhang .  -mat_mumps_icntl_3 <0>: ICNTL(3): output stream for global information, collected on the host (None)
22834e34a73bSHong Zhang .  -mat_mumps_icntl_4 <0>: ICNTL(4): level of printing (0 to 4) (None)
22844e34a73bSHong Zhang .  -mat_mumps_icntl_6 <7>: ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7) (None)
22854e34a73bSHong 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)
22864e34a73bSHong Zhang .  -mat_mumps_icntl_8 <77>: ICNTL(8): scaling strategy (-2 to 8 or 77) (None)
22874e34a73bSHong Zhang .  -mat_mumps_icntl_10 <0>: ICNTL(10): max num of refinements (None)
22884e34a73bSHong Zhang .  -mat_mumps_icntl_11 <0>: ICNTL(11): statistics related to an error analysis (via -ksp_view) (None)
22894e34a73bSHong Zhang .  -mat_mumps_icntl_12 <1>: ICNTL(12): an ordering strategy for symmetric matrices (0 to 3) (None)
22904e34a73bSHong Zhang .  -mat_mumps_icntl_13 <0>: ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting (None)
22914e34a73bSHong Zhang .  -mat_mumps_icntl_14 <20>: ICNTL(14): percentage increase in the estimated working space (None)
22924e34a73bSHong Zhang .  -mat_mumps_icntl_19 <0>: ICNTL(19): computes the Schur complement (None)
22934e34a73bSHong Zhang .  -mat_mumps_icntl_22 <0>: ICNTL(22): in-core/out-of-core factorization and solve (0 or 1) (None)
22944e34a73bSHong Zhang .  -mat_mumps_icntl_23 <0>: ICNTL(23): max size of the working memory (MB) that can allocate per processor (None)
22954e34a73bSHong Zhang .  -mat_mumps_icntl_24 <0>: ICNTL(24): detection of null pivot rows (0 or 1) (None)
22964e34a73bSHong Zhang .  -mat_mumps_icntl_25 <0>: ICNTL(25): compute a solution of a deficient matrix and a null space basis (None)
22974e34a73bSHong Zhang .  -mat_mumps_icntl_26 <0>: ICNTL(26): drives the solution phase if a Schur complement matrix (None)
22984e34a73bSHong 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)
22994e34a73bSHong Zhang .  -mat_mumps_icntl_29 <0>: ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis (None)
23004e34a73bSHong Zhang .  -mat_mumps_icntl_30 <0>: ICNTL(30): compute user-specified set of entries in inv(A) (None)
23014e34a73bSHong Zhang .  -mat_mumps_icntl_31 <0>: ICNTL(31): indicates which factors may be discarded during factorization (None)
23024e34a73bSHong Zhang .  -mat_mumps_icntl_33 <0>: ICNTL(33): compute determinant (None)
23034e34a73bSHong Zhang .  -mat_mumps_cntl_1 <0.01>: CNTL(1): relative pivoting threshold (None)
23044e34a73bSHong Zhang .  -mat_mumps_cntl_2 <1.49012e-08>: CNTL(2): stopping criterion of refinement (None)
23054e34a73bSHong Zhang .  -mat_mumps_cntl_3 <0>: CNTL(3): absolute pivoting threshold (None)
23064e34a73bSHong Zhang .  -mat_mumps_cntl_4 <-1>: CNTL(4): value for static pivoting (None)
23074e34a73bSHong Zhang -  -mat_mumps_cntl_5 <0>: CNTL(5): fixation for null pivots (None)
230824b6179bSKris Buschelman 
230924b6179bSKris Buschelman   Level: beginner
231024b6179bSKris Buschelman 
231141c8de11SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage
231241c8de11SBarry Smith 
231324b6179bSKris Buschelman M*/
231424b6179bSKris Buschelman 
231535bd34faSBarry Smith #undef __FUNCT__
231635bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_mumps"
2317f7a08781SBarry Smith static PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type)
231835bd34faSBarry Smith {
231935bd34faSBarry Smith   PetscFunctionBegin;
23202692d6eeSBarry Smith   *type = MATSOLVERMUMPS;
232135bd34faSBarry Smith   PetscFunctionReturn(0);
232235bd34faSBarry Smith }
232335bd34faSBarry Smith 
2324bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI AIJ matrices */
23252877fffaSHong Zhang #undef __FUNCT__
2326bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_aij_mumps"
23278cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
23282877fffaSHong Zhang {
23292877fffaSHong Zhang   Mat            B;
23302877fffaSHong Zhang   PetscErrorCode ierr;
23312877fffaSHong Zhang   Mat_MUMPS      *mumps;
2332ace3abfcSBarry Smith   PetscBool      isSeqAIJ;
23332877fffaSHong Zhang 
23342877fffaSHong Zhang   PetscFunctionBegin;
23352877fffaSHong Zhang   /* Create the factorization matrix */
2336251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr);
2337ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
23382877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
23392877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
2340bccb9932SShri Abhyankar   if (isSeqAIJ) {
23410298fd71SBarry Smith     ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr);
2342bccb9932SShri Abhyankar   } else {
23430298fd71SBarry Smith     ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr);
2344bccb9932SShri Abhyankar   }
23452877fffaSHong Zhang 
2346b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
23472205254eSKarl Rupp 
23482877fffaSHong Zhang   B->ops->view        = MatView_MUMPS;
234935bd34faSBarry Smith   B->ops->getinfo     = MatGetInfo_MUMPS;
235020be8e61SHong Zhang   B->ops->getdiagonal = MatGetDiagonal_MUMPS;
23512205254eSKarl Rupp 
2352bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
23535a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr);
23545a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorInvertSchurComplement_C",MatFactorInvertSchurComplement_MUMPS);CHKERRQ(ierr);
23555a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr);
23565a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSchurComplement_C",MatFactorGetSchurComplement_MUMPS);CHKERRQ(ierr);
23575a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplement_C",MatFactorSolveSchurComplement_MUMPS);CHKERRQ(ierr);
23585a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplementTranspose_C",MatFactorSolveSchurComplementTranspose_MUMPS);CHKERRQ(ierr);
2359bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
2360bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
2361bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
2362bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
2363ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
2364ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
2365ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
2366ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
23676444a565SStefano Zampini 
2368450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
2369450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
2370d5f3da31SBarry Smith     B->factortype            = MAT_FACTOR_LU;
2371bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
2372bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
2373746480a1SHong Zhang     mumps->sym = 0;
2374dcd589f8SShri Abhyankar   } else {
237567877ebaSShri Abhyankar     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2376450b117fSShri Abhyankar     B->factortype                  = MAT_FACTOR_CHOLESKY;
2377bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
2378bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
237959ac8732SStefano Zampini #if defined(PETSC_USE_COMPLEX)
238059ac8732SStefano Zampini     mumps->sym = 2;
238159ac8732SStefano Zampini #else
23826fdc2a6dSBarry Smith     if (A->spd_set && A->spd) mumps->sym = 1;
23836fdc2a6dSBarry Smith     else                      mumps->sym = 2;
238459ac8732SStefano Zampini #endif
2385450b117fSShri Abhyankar   }
23862877fffaSHong Zhang 
23872877fffaSHong Zhang   mumps->isAIJ    = PETSC_TRUE;
2388bf0cc555SLisandro Dalcin   mumps->Destroy  = B->ops->destroy;
23892877fffaSHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
23902877fffaSHong Zhang   B->spptr        = (void*)mumps;
23912205254eSKarl Rupp 
2392f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
2393746480a1SHong Zhang 
23942877fffaSHong Zhang   *F = B;
23952877fffaSHong Zhang   PetscFunctionReturn(0);
23962877fffaSHong Zhang }
23972877fffaSHong Zhang 
2398bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI SBAIJ matrices */
23992877fffaSHong Zhang #undef __FUNCT__
2400bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_sbaij_mumps"
24018cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
24022877fffaSHong Zhang {
24032877fffaSHong Zhang   Mat            B;
24042877fffaSHong Zhang   PetscErrorCode ierr;
24052877fffaSHong Zhang   Mat_MUMPS      *mumps;
2406ace3abfcSBarry Smith   PetscBool      isSeqSBAIJ;
24072877fffaSHong Zhang 
24082877fffaSHong Zhang   PetscFunctionBegin;
2409ce94432eSBarry Smith   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
2410ce94432eSBarry 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");
2411251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
24122877fffaSHong Zhang   /* Create the factorization matrix */
2413ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
24142877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
24152877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
2416b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
2417bccb9932SShri Abhyankar   if (isSeqSBAIJ) {
24180298fd71SBarry Smith     ierr = MatSeqSBAIJSetPreallocation(B,1,0,NULL);CHKERRQ(ierr);
24192205254eSKarl Rupp 
242016ebf90aSShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
2421dcd589f8SShri Abhyankar   } else {
24220298fd71SBarry Smith     ierr = MatMPISBAIJSetPreallocation(B,1,0,NULL,0,NULL);CHKERRQ(ierr);
24232205254eSKarl Rupp 
2424bccb9932SShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
2425bccb9932SShri Abhyankar   }
2426bccb9932SShri Abhyankar 
242767877ebaSShri Abhyankar   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2428bccb9932SShri Abhyankar   B->ops->view                   = MatView_MUMPS;
242920be8e61SHong Zhang   B->ops->getdiagonal            = MatGetDiagonal_MUMPS;
24302205254eSKarl Rupp 
2431bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
24325a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr);
24335a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorInvertSchurComplement_C",MatFactorInvertSchurComplement_MUMPS);CHKERRQ(ierr);
24345a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr);
24355a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSchurComplement_C",MatFactorGetSchurComplement_MUMPS);CHKERRQ(ierr);
24365a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplement_C",MatFactorSolveSchurComplement_MUMPS);CHKERRQ(ierr);
24375a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplementTranspose_C",MatFactorSolveSchurComplementTranspose_MUMPS);CHKERRQ(ierr);
2438b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
2439b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
2440b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
2441b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
2442ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
2443ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
2444ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
2445ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
24462205254eSKarl Rupp 
2447f4762488SHong Zhang   B->factortype = MAT_FACTOR_CHOLESKY;
244859ac8732SStefano Zampini #if defined(PETSC_USE_COMPLEX)
244959ac8732SStefano Zampini   mumps->sym = 2;
245059ac8732SStefano Zampini #else
24516fdc2a6dSBarry Smith   if (A->spd_set && A->spd) mumps->sym = 1;
24526fdc2a6dSBarry Smith   else                      mumps->sym = 2;
245359ac8732SStefano Zampini #endif
2454a214ac2aSShri Abhyankar 
2455bccb9932SShri Abhyankar   mumps->isAIJ    = PETSC_FALSE;
2456bf0cc555SLisandro Dalcin   mumps->Destroy  = B->ops->destroy;
2457f3c0ef26SHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
24582877fffaSHong Zhang   B->spptr        = (void*)mumps;
24592205254eSKarl Rupp 
2460f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
2461746480a1SHong Zhang 
24622877fffaSHong Zhang   *F = B;
24632877fffaSHong Zhang   PetscFunctionReturn(0);
24642877fffaSHong Zhang }
246597969023SHong Zhang 
2466450b117fSShri Abhyankar #undef __FUNCT__
2467bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_baij_mumps"
24688cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
246967877ebaSShri Abhyankar {
247067877ebaSShri Abhyankar   Mat            B;
247167877ebaSShri Abhyankar   PetscErrorCode ierr;
247267877ebaSShri Abhyankar   Mat_MUMPS      *mumps;
2473ace3abfcSBarry Smith   PetscBool      isSeqBAIJ;
247467877ebaSShri Abhyankar 
247567877ebaSShri Abhyankar   PetscFunctionBegin;
247667877ebaSShri Abhyankar   /* Create the factorization matrix */
2477251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr);
2478ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
247967877ebaSShri Abhyankar   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
248067877ebaSShri Abhyankar   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
2481bccb9932SShri Abhyankar   if (isSeqBAIJ) {
24820298fd71SBarry Smith     ierr = MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,NULL);CHKERRQ(ierr);
2483bccb9932SShri Abhyankar   } else {
24840298fd71SBarry Smith     ierr = MatMPIBAIJSetPreallocation(B,A->rmap->bs,0,NULL,0,NULL);CHKERRQ(ierr);
2485bccb9932SShri Abhyankar   }
2486450b117fSShri Abhyankar 
2487b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
2488450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
2489450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
2490450b117fSShri Abhyankar     B->factortype            = MAT_FACTOR_LU;
2491bccb9932SShri Abhyankar     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
2492bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
2493746480a1SHong Zhang     mumps->sym = 0;
2494f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n");
2495bccb9932SShri Abhyankar 
2496450b117fSShri Abhyankar   B->ops->view        = MatView_MUMPS;
249720be8e61SHong Zhang   B->ops->getdiagonal = MatGetDiagonal_MUMPS;
24982205254eSKarl Rupp 
2499bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
25005a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr);
25015a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorInvertSchurComplement_C",MatFactorInvertSchurComplement_MUMPS);CHKERRQ(ierr);
25025a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr);
25035a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSchurComplement_C",MatFactorGetSchurComplement_MUMPS);CHKERRQ(ierr);
25045a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplement_C",MatFactorSolveSchurComplement_MUMPS);CHKERRQ(ierr);
25055a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplementTranspose_C",MatFactorSolveSchurComplementTranspose_MUMPS);CHKERRQ(ierr);
2506bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
2507bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
2508bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
2509bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
2510ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
2511ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
2512ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
2513ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
2514450b117fSShri Abhyankar 
2515450b117fSShri Abhyankar   mumps->isAIJ    = PETSC_TRUE;
2516bf0cc555SLisandro Dalcin   mumps->Destroy  = B->ops->destroy;
2517450b117fSShri Abhyankar   B->ops->destroy = MatDestroy_MUMPS;
2518450b117fSShri Abhyankar   B->spptr        = (void*)mumps;
25192205254eSKarl Rupp 
2520f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
2521746480a1SHong Zhang 
2522450b117fSShri Abhyankar   *F = B;
2523450b117fSShri Abhyankar   PetscFunctionReturn(0);
2524450b117fSShri Abhyankar }
252542c9c57cSBarry Smith 
252642c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mumps(Mat,MatFactorType,Mat*);
252742c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_baij_mumps(Mat,MatFactorType,Mat*);
252842c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat,MatFactorType,Mat*);
252942c9c57cSBarry Smith 
253042c9c57cSBarry Smith #undef __FUNCT__
253142c9c57cSBarry Smith #define __FUNCT__ "MatSolverPackageRegister_MUMPS"
253229b38603SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_MUMPS(void)
253342c9c57cSBarry Smith {
253442c9c57cSBarry Smith   PetscErrorCode ierr;
253542c9c57cSBarry Smith 
253642c9c57cSBarry Smith   PetscFunctionBegin;
253742c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr);
253842c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr);
253942c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr);
254042c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr);
254142c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPISBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr);
254242c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr);
254342c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr);
254442c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr);
254542c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr);
254642c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr);
254742c9c57cSBarry Smith   PetscFunctionReturn(0);
254842c9c57cSBarry Smith }
254942c9c57cSBarry Smith 
2550