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 PetscInt sizeredrhs = mumps->id.nrhs*mumps->id.size_schur; 348b5fa320bSStefano Zampini /* allocate MUMPS internal array to store reduced right-hand sides */ 349b5fa320bSStefano Zampini if (!mumps->id.redrhs || sizeredrhs > mumps->sizeredrhs) { 350b5fa320bSStefano Zampini ierr = PetscFree(mumps->id.redrhs);CHKERRQ(ierr); 351b5fa320bSStefano Zampini mumps->id.lredrhs = mumps->id.size_schur; 352b5fa320bSStefano Zampini ierr = PetscMalloc1(mumps->id.nrhs*mumps->id.lredrhs,&mumps->id.redrhs);CHKERRQ(ierr); 353b5fa320bSStefano Zampini mumps->sizeredrhs = mumps->id.nrhs*mumps->id.lredrhs; 354b5fa320bSStefano Zampini } 355b5fa320bSStefano Zampini mumps->id.ICNTL(26) = 1; /* condensation phase */ 356b5fa320bSStefano Zampini } else { /* prepare for the expansion step */ 357b8f61ee1SStefano Zampini /* solve Schur complement (this has to be done by the MUMPS user, so basically us) */ 358e807eca7SStefano Zampini ierr = MatMumpsSolveSchur_Private(mumps,PETSC_TRUE);CHKERRQ(ierr); 359b5fa320bSStefano Zampini mumps->id.ICNTL(26) = 2; /* expansion phase */ 360b5fa320bSStefano Zampini PetscMUMPS_c(&mumps->id); 361b5fa320bSStefano 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)); 362b5fa320bSStefano Zampini /* restore defaults */ 363b5fa320bSStefano Zampini mumps->id.ICNTL(26) = -1; 364b5fa320bSStefano Zampini } 365b5fa320bSStefano Zampini PetscFunctionReturn(0); 366b5fa320bSStefano Zampini } 367b5fa320bSStefano Zampini 368397b6df1SKris Buschelman /* 369d341cd04SHong Zhang MatConvertToTriples_A_B - convert Petsc matrix to triples: row[nz], col[nz], val[nz] 370d341cd04SHong Zhang 371397b6df1SKris Buschelman input: 37267877ebaSShri Abhyankar A - matrix in aij,baij or sbaij (bs=1) format 373397b6df1SKris Buschelman shift - 0: C style output triple; 1: Fortran style output triple. 374bccb9932SShri Abhyankar reuse - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple 375bccb9932SShri Abhyankar MAT_REUSE_MATRIX: only the values in v array are updated 376397b6df1SKris Buschelman output: 377397b6df1SKris Buschelman nnz - dim of r, c, and v (number of local nonzero entries of A) 378397b6df1SKris Buschelman r, c, v - row and col index, matrix values (matrix triples) 379eb9baa12SBarry Smith 380eb9baa12SBarry Smith The returned values r, c, and sometimes v are obtained in a single PetscMalloc(). Then in MatDestroy_MUMPS() it is 381eb9baa12SBarry Smith freed with PetscFree((mumps->irn); This is not ideal code, the fact that v is ONLY sometimes part of mumps->irn means 382eb9baa12SBarry Smith that the PetscMalloc() cannot easily be replaced with a PetscMalloc3(). 383eb9baa12SBarry Smith 384397b6df1SKris Buschelman */ 38516ebf90aSShri Abhyankar 38616ebf90aSShri Abhyankar #undef __FUNCT__ 38716ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqaij" 388bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 389b24902e0SBarry Smith { 390185f6596SHong Zhang const PetscInt *ai,*aj,*ajj,M=A->rmap->n; 39167877ebaSShri Abhyankar PetscInt nz,rnz,i,j; 392dfbe8321SBarry Smith PetscErrorCode ierr; 393c1490034SHong Zhang PetscInt *row,*col; 39416ebf90aSShri Abhyankar Mat_SeqAIJ *aa=(Mat_SeqAIJ*)A->data; 395397b6df1SKris Buschelman 396397b6df1SKris Buschelman PetscFunctionBegin; 39716ebf90aSShri Abhyankar *v=aa->a; 398bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 3992205254eSKarl Rupp nz = aa->nz; 4002205254eSKarl Rupp ai = aa->i; 4012205254eSKarl Rupp aj = aa->j; 40216ebf90aSShri Abhyankar *nnz = nz; 403785e854fSJed Brown ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr); 404185f6596SHong Zhang col = row + nz; 405185f6596SHong Zhang 40616ebf90aSShri Abhyankar nz = 0; 40716ebf90aSShri Abhyankar for (i=0; i<M; i++) { 40816ebf90aSShri Abhyankar rnz = ai[i+1] - ai[i]; 40967877ebaSShri Abhyankar ajj = aj + ai[i]; 41067877ebaSShri Abhyankar for (j=0; j<rnz; j++) { 41167877ebaSShri Abhyankar row[nz] = i+shift; col[nz++] = ajj[j] + shift; 41216ebf90aSShri Abhyankar } 41316ebf90aSShri Abhyankar } 41416ebf90aSShri Abhyankar *r = row; *c = col; 41516ebf90aSShri Abhyankar } 41616ebf90aSShri Abhyankar PetscFunctionReturn(0); 41716ebf90aSShri Abhyankar } 418397b6df1SKris Buschelman 41916ebf90aSShri Abhyankar #undef __FUNCT__ 42067877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqbaij_seqaij" 421bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 42267877ebaSShri Abhyankar { 42367877ebaSShri Abhyankar Mat_SeqBAIJ *aa=(Mat_SeqBAIJ*)A->data; 42433d57670SJed Brown const PetscInt *ai,*aj,*ajj,bs2 = aa->bs2; 42533d57670SJed Brown PetscInt bs,M,nz,idx=0,rnz,i,j,k,m; 42667877ebaSShri Abhyankar PetscErrorCode ierr; 42767877ebaSShri Abhyankar PetscInt *row,*col; 42867877ebaSShri Abhyankar 42967877ebaSShri Abhyankar PetscFunctionBegin; 43033d57670SJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 43133d57670SJed Brown M = A->rmap->N/bs; 432cf3759fdSShri Abhyankar *v = aa->a; 433bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 434cf3759fdSShri Abhyankar ai = aa->i; aj = aa->j; 43567877ebaSShri Abhyankar nz = bs2*aa->nz; 43667877ebaSShri Abhyankar *nnz = nz; 437785e854fSJed Brown ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr); 438185f6596SHong Zhang col = row + nz; 439185f6596SHong Zhang 44067877ebaSShri Abhyankar for (i=0; i<M; i++) { 44167877ebaSShri Abhyankar ajj = aj + ai[i]; 44267877ebaSShri Abhyankar rnz = ai[i+1] - ai[i]; 44367877ebaSShri Abhyankar for (k=0; k<rnz; k++) { 44467877ebaSShri Abhyankar for (j=0; j<bs; j++) { 44567877ebaSShri Abhyankar for (m=0; m<bs; m++) { 44667877ebaSShri Abhyankar row[idx] = i*bs + m + shift; 447cf3759fdSShri Abhyankar col[idx++] = bs*(ajj[k]) + j + shift; 44867877ebaSShri Abhyankar } 44967877ebaSShri Abhyankar } 45067877ebaSShri Abhyankar } 45167877ebaSShri Abhyankar } 452cf3759fdSShri Abhyankar *r = row; *c = col; 45367877ebaSShri Abhyankar } 45467877ebaSShri Abhyankar PetscFunctionReturn(0); 45567877ebaSShri Abhyankar } 45667877ebaSShri Abhyankar 45767877ebaSShri Abhyankar #undef __FUNCT__ 45816ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqsbaij_seqsbaij" 459bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 46016ebf90aSShri Abhyankar { 46167877ebaSShri Abhyankar const PetscInt *ai, *aj,*ajj,M=A->rmap->n; 46267877ebaSShri Abhyankar PetscInt nz,rnz,i,j; 46316ebf90aSShri Abhyankar PetscErrorCode ierr; 46416ebf90aSShri Abhyankar PetscInt *row,*col; 46516ebf90aSShri Abhyankar Mat_SeqSBAIJ *aa=(Mat_SeqSBAIJ*)A->data; 46616ebf90aSShri Abhyankar 46716ebf90aSShri Abhyankar PetscFunctionBegin; 468882afa5aSHong Zhang *v = aa->a; 469bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 4702205254eSKarl Rupp nz = aa->nz; 4712205254eSKarl Rupp ai = aa->i; 4722205254eSKarl Rupp aj = aa->j; 4732205254eSKarl Rupp *v = aa->a; 47416ebf90aSShri Abhyankar *nnz = nz; 475785e854fSJed Brown ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr); 476185f6596SHong Zhang col = row + nz; 477185f6596SHong Zhang 47816ebf90aSShri Abhyankar nz = 0; 47916ebf90aSShri Abhyankar for (i=0; i<M; i++) { 48016ebf90aSShri Abhyankar rnz = ai[i+1] - ai[i]; 48167877ebaSShri Abhyankar ajj = aj + ai[i]; 48267877ebaSShri Abhyankar for (j=0; j<rnz; j++) { 48367877ebaSShri Abhyankar row[nz] = i+shift; col[nz++] = ajj[j] + shift; 48416ebf90aSShri Abhyankar } 48516ebf90aSShri Abhyankar } 48616ebf90aSShri Abhyankar *r = row; *c = col; 48716ebf90aSShri Abhyankar } 48816ebf90aSShri Abhyankar PetscFunctionReturn(0); 48916ebf90aSShri Abhyankar } 49016ebf90aSShri Abhyankar 49116ebf90aSShri Abhyankar #undef __FUNCT__ 49216ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqsbaij" 493bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 49416ebf90aSShri Abhyankar { 49567877ebaSShri Abhyankar const PetscInt *ai,*aj,*ajj,*adiag,M=A->rmap->n; 49667877ebaSShri Abhyankar PetscInt nz,rnz,i,j; 49767877ebaSShri Abhyankar const PetscScalar *av,*v1; 49816ebf90aSShri Abhyankar PetscScalar *val; 49916ebf90aSShri Abhyankar PetscErrorCode ierr; 50016ebf90aSShri Abhyankar PetscInt *row,*col; 501829b1710SHong Zhang Mat_SeqAIJ *aa=(Mat_SeqAIJ*)A->data; 50216ebf90aSShri Abhyankar 50316ebf90aSShri Abhyankar PetscFunctionBegin; 50416ebf90aSShri Abhyankar ai =aa->i; aj=aa->j;av=aa->a; 50516ebf90aSShri Abhyankar adiag=aa->diag; 506bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 507829b1710SHong Zhang /* count nz in the uppper triangular part of A */ 508829b1710SHong Zhang nz = 0; 509829b1710SHong Zhang for (i=0; i<M; i++) nz += ai[i+1] - adiag[i]; 51016ebf90aSShri Abhyankar *nnz = nz; 511829b1710SHong Zhang 512185f6596SHong Zhang ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 513185f6596SHong Zhang col = row + nz; 514185f6596SHong Zhang val = (PetscScalar*)(col + nz); 515185f6596SHong Zhang 51616ebf90aSShri Abhyankar nz = 0; 51716ebf90aSShri Abhyankar for (i=0; i<M; i++) { 51816ebf90aSShri Abhyankar rnz = ai[i+1] - adiag[i]; 51967877ebaSShri Abhyankar ajj = aj + adiag[i]; 520cf3759fdSShri Abhyankar v1 = av + adiag[i]; 52167877ebaSShri Abhyankar for (j=0; j<rnz; j++) { 52267877ebaSShri Abhyankar row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j]; 52316ebf90aSShri Abhyankar } 52416ebf90aSShri Abhyankar } 52516ebf90aSShri Abhyankar *r = row; *c = col; *v = val; 526397b6df1SKris Buschelman } else { 52716ebf90aSShri Abhyankar nz = 0; val = *v; 52816ebf90aSShri Abhyankar for (i=0; i <M; i++) { 52916ebf90aSShri Abhyankar rnz = ai[i+1] - adiag[i]; 53067877ebaSShri Abhyankar ajj = aj + adiag[i]; 53167877ebaSShri Abhyankar v1 = av + adiag[i]; 53267877ebaSShri Abhyankar for (j=0; j<rnz; j++) { 53367877ebaSShri Abhyankar val[nz++] = v1[j]; 53416ebf90aSShri Abhyankar } 53516ebf90aSShri Abhyankar } 53616ebf90aSShri Abhyankar } 53716ebf90aSShri Abhyankar PetscFunctionReturn(0); 53816ebf90aSShri Abhyankar } 53916ebf90aSShri Abhyankar 54016ebf90aSShri Abhyankar #undef __FUNCT__ 54116ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpisbaij_mpisbaij" 542bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 54316ebf90aSShri Abhyankar { 54416ebf90aSShri Abhyankar const PetscInt *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj; 54516ebf90aSShri Abhyankar PetscErrorCode ierr; 54616ebf90aSShri Abhyankar PetscInt rstart,nz,i,j,jj,irow,countA,countB; 54716ebf90aSShri Abhyankar PetscInt *row,*col; 54816ebf90aSShri Abhyankar const PetscScalar *av, *bv,*v1,*v2; 54916ebf90aSShri Abhyankar PetscScalar *val; 550397b6df1SKris Buschelman Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)A->data; 551397b6df1SKris Buschelman Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)(mat->A)->data; 552397b6df1SKris Buschelman Mat_SeqBAIJ *bb = (Mat_SeqBAIJ*)(mat->B)->data; 55316ebf90aSShri Abhyankar 55416ebf90aSShri Abhyankar PetscFunctionBegin; 555d0f46423SBarry Smith ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart; 556397b6df1SKris Buschelman av=aa->a; bv=bb->a; 557397b6df1SKris Buschelman 5582205254eSKarl Rupp garray = mat->garray; 5592205254eSKarl Rupp 560bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 56116ebf90aSShri Abhyankar nz = aa->nz + bb->nz; 56216ebf90aSShri Abhyankar *nnz = nz; 563185f6596SHong Zhang ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 564185f6596SHong Zhang col = row + nz; 565185f6596SHong Zhang val = (PetscScalar*)(col + nz); 566185f6596SHong Zhang 567397b6df1SKris Buschelman *r = row; *c = col; *v = val; 568397b6df1SKris Buschelman } else { 569397b6df1SKris Buschelman row = *r; col = *c; val = *v; 570397b6df1SKris Buschelman } 571397b6df1SKris Buschelman 572028e57e8SHong Zhang jj = 0; irow = rstart; 573397b6df1SKris Buschelman for (i=0; i<m; i++) { 574397b6df1SKris Buschelman ajj = aj + ai[i]; /* ptr to the beginning of this row */ 575397b6df1SKris Buschelman countA = ai[i+1] - ai[i]; 576397b6df1SKris Buschelman countB = bi[i+1] - bi[i]; 577397b6df1SKris Buschelman bjj = bj + bi[i]; 57816ebf90aSShri Abhyankar v1 = av + ai[i]; 57916ebf90aSShri Abhyankar v2 = bv + bi[i]; 580397b6df1SKris Buschelman 581397b6df1SKris Buschelman /* A-part */ 582397b6df1SKris Buschelman for (j=0; j<countA; j++) { 583bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 584397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 585397b6df1SKris Buschelman } 58616ebf90aSShri Abhyankar val[jj++] = v1[j]; 587397b6df1SKris Buschelman } 58816ebf90aSShri Abhyankar 58916ebf90aSShri Abhyankar /* B-part */ 59016ebf90aSShri Abhyankar for (j=0; j < countB; j++) { 591bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 592397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 593397b6df1SKris Buschelman } 59416ebf90aSShri Abhyankar val[jj++] = v2[j]; 59516ebf90aSShri Abhyankar } 59616ebf90aSShri Abhyankar irow++; 59716ebf90aSShri Abhyankar } 59816ebf90aSShri Abhyankar PetscFunctionReturn(0); 59916ebf90aSShri Abhyankar } 60016ebf90aSShri Abhyankar 60116ebf90aSShri Abhyankar #undef __FUNCT__ 60216ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpiaij" 603bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 60416ebf90aSShri Abhyankar { 60516ebf90aSShri Abhyankar const PetscInt *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj; 60616ebf90aSShri Abhyankar PetscErrorCode ierr; 60716ebf90aSShri Abhyankar PetscInt rstart,nz,i,j,jj,irow,countA,countB; 60816ebf90aSShri Abhyankar PetscInt *row,*col; 60916ebf90aSShri Abhyankar const PetscScalar *av, *bv,*v1,*v2; 61016ebf90aSShri Abhyankar PetscScalar *val; 61116ebf90aSShri Abhyankar Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data; 61216ebf90aSShri Abhyankar Mat_SeqAIJ *aa = (Mat_SeqAIJ*)(mat->A)->data; 61316ebf90aSShri Abhyankar Mat_SeqAIJ *bb = (Mat_SeqAIJ*)(mat->B)->data; 61416ebf90aSShri Abhyankar 61516ebf90aSShri Abhyankar PetscFunctionBegin; 61616ebf90aSShri Abhyankar ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart; 61716ebf90aSShri Abhyankar av=aa->a; bv=bb->a; 61816ebf90aSShri Abhyankar 6192205254eSKarl Rupp garray = mat->garray; 6202205254eSKarl Rupp 621bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 62216ebf90aSShri Abhyankar nz = aa->nz + bb->nz; 62316ebf90aSShri Abhyankar *nnz = nz; 624185f6596SHong Zhang ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 625185f6596SHong Zhang col = row + nz; 626185f6596SHong Zhang val = (PetscScalar*)(col + nz); 627185f6596SHong Zhang 62816ebf90aSShri Abhyankar *r = row; *c = col; *v = val; 62916ebf90aSShri Abhyankar } else { 63016ebf90aSShri Abhyankar row = *r; col = *c; val = *v; 63116ebf90aSShri Abhyankar } 63216ebf90aSShri Abhyankar 63316ebf90aSShri Abhyankar jj = 0; irow = rstart; 63416ebf90aSShri Abhyankar for (i=0; i<m; i++) { 63516ebf90aSShri Abhyankar ajj = aj + ai[i]; /* ptr to the beginning of this row */ 63616ebf90aSShri Abhyankar countA = ai[i+1] - ai[i]; 63716ebf90aSShri Abhyankar countB = bi[i+1] - bi[i]; 63816ebf90aSShri Abhyankar bjj = bj + bi[i]; 63916ebf90aSShri Abhyankar v1 = av + ai[i]; 64016ebf90aSShri Abhyankar v2 = bv + bi[i]; 64116ebf90aSShri Abhyankar 64216ebf90aSShri Abhyankar /* A-part */ 64316ebf90aSShri Abhyankar for (j=0; j<countA; j++) { 644bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 64516ebf90aSShri Abhyankar row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 64616ebf90aSShri Abhyankar } 64716ebf90aSShri Abhyankar val[jj++] = v1[j]; 64816ebf90aSShri Abhyankar } 64916ebf90aSShri Abhyankar 65016ebf90aSShri Abhyankar /* B-part */ 65116ebf90aSShri Abhyankar for (j=0; j < countB; j++) { 652bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 65316ebf90aSShri Abhyankar row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 65416ebf90aSShri Abhyankar } 65516ebf90aSShri Abhyankar val[jj++] = v2[j]; 65616ebf90aSShri Abhyankar } 65716ebf90aSShri Abhyankar irow++; 65816ebf90aSShri Abhyankar } 65916ebf90aSShri Abhyankar PetscFunctionReturn(0); 66016ebf90aSShri Abhyankar } 66116ebf90aSShri Abhyankar 66216ebf90aSShri Abhyankar #undef __FUNCT__ 66367877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpibaij_mpiaij" 664bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 66567877ebaSShri Abhyankar { 66667877ebaSShri Abhyankar Mat_MPIBAIJ *mat = (Mat_MPIBAIJ*)A->data; 66767877ebaSShri Abhyankar Mat_SeqBAIJ *aa = (Mat_SeqBAIJ*)(mat->A)->data; 66867877ebaSShri Abhyankar Mat_SeqBAIJ *bb = (Mat_SeqBAIJ*)(mat->B)->data; 66967877ebaSShri Abhyankar const PetscInt *ai = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj; 670d985c460SShri Abhyankar const PetscInt *garray = mat->garray,mbs=mat->mbs,rstart=A->rmap->rstart; 67133d57670SJed Brown const PetscInt bs2=mat->bs2; 67267877ebaSShri Abhyankar PetscErrorCode ierr; 67333d57670SJed Brown PetscInt bs,nz,i,j,k,n,jj,irow,countA,countB,idx; 67467877ebaSShri Abhyankar PetscInt *row,*col; 67567877ebaSShri Abhyankar const PetscScalar *av=aa->a, *bv=bb->a,*v1,*v2; 67667877ebaSShri Abhyankar PetscScalar *val; 67767877ebaSShri Abhyankar 67867877ebaSShri Abhyankar PetscFunctionBegin; 67933d57670SJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 680bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 68167877ebaSShri Abhyankar nz = bs2*(aa->nz + bb->nz); 68267877ebaSShri Abhyankar *nnz = nz; 683185f6596SHong Zhang ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 684185f6596SHong Zhang col = row + nz; 685185f6596SHong Zhang val = (PetscScalar*)(col + nz); 686185f6596SHong Zhang 68767877ebaSShri Abhyankar *r = row; *c = col; *v = val; 68867877ebaSShri Abhyankar } else { 68967877ebaSShri Abhyankar row = *r; col = *c; val = *v; 69067877ebaSShri Abhyankar } 69167877ebaSShri Abhyankar 692d985c460SShri Abhyankar jj = 0; irow = rstart; 69367877ebaSShri Abhyankar for (i=0; i<mbs; i++) { 69467877ebaSShri Abhyankar countA = ai[i+1] - ai[i]; 69567877ebaSShri Abhyankar countB = bi[i+1] - bi[i]; 69667877ebaSShri Abhyankar ajj = aj + ai[i]; 69767877ebaSShri Abhyankar bjj = bj + bi[i]; 69867877ebaSShri Abhyankar v1 = av + bs2*ai[i]; 69967877ebaSShri Abhyankar v2 = bv + bs2*bi[i]; 70067877ebaSShri Abhyankar 70167877ebaSShri Abhyankar idx = 0; 70267877ebaSShri Abhyankar /* A-part */ 70367877ebaSShri Abhyankar for (k=0; k<countA; k++) { 70467877ebaSShri Abhyankar for (j=0; j<bs; j++) { 70567877ebaSShri Abhyankar for (n=0; n<bs; n++) { 706bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 707d985c460SShri Abhyankar row[jj] = irow + n + shift; 708d985c460SShri Abhyankar col[jj] = rstart + bs*ajj[k] + j + shift; 70967877ebaSShri Abhyankar } 71067877ebaSShri Abhyankar val[jj++] = v1[idx++]; 71167877ebaSShri Abhyankar } 71267877ebaSShri Abhyankar } 71367877ebaSShri Abhyankar } 71467877ebaSShri Abhyankar 71567877ebaSShri Abhyankar idx = 0; 71667877ebaSShri Abhyankar /* B-part */ 71767877ebaSShri Abhyankar for (k=0; k<countB; k++) { 71867877ebaSShri Abhyankar for (j=0; j<bs; j++) { 71967877ebaSShri Abhyankar for (n=0; n<bs; n++) { 720bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 721d985c460SShri Abhyankar row[jj] = irow + n + shift; 722d985c460SShri Abhyankar col[jj] = bs*garray[bjj[k]] + j + shift; 72367877ebaSShri Abhyankar } 724d985c460SShri Abhyankar val[jj++] = v2[idx++]; 72567877ebaSShri Abhyankar } 72667877ebaSShri Abhyankar } 72767877ebaSShri Abhyankar } 728d985c460SShri Abhyankar irow += bs; 72967877ebaSShri Abhyankar } 73067877ebaSShri Abhyankar PetscFunctionReturn(0); 73167877ebaSShri Abhyankar } 73267877ebaSShri Abhyankar 73367877ebaSShri Abhyankar #undef __FUNCT__ 73416ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpisbaij" 735bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 73616ebf90aSShri Abhyankar { 73716ebf90aSShri Abhyankar const PetscInt *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj; 73816ebf90aSShri Abhyankar PetscErrorCode ierr; 739e0bace9bSHong Zhang PetscInt rstart,nz,nza,nzb,i,j,jj,irow,countA,countB; 74016ebf90aSShri Abhyankar PetscInt *row,*col; 74116ebf90aSShri Abhyankar const PetscScalar *av, *bv,*v1,*v2; 74216ebf90aSShri Abhyankar PetscScalar *val; 74316ebf90aSShri Abhyankar Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data; 74416ebf90aSShri Abhyankar Mat_SeqAIJ *aa =(Mat_SeqAIJ*)(mat->A)->data; 74516ebf90aSShri Abhyankar Mat_SeqAIJ *bb =(Mat_SeqAIJ*)(mat->B)->data; 74616ebf90aSShri Abhyankar 74716ebf90aSShri Abhyankar PetscFunctionBegin; 74816ebf90aSShri Abhyankar ai=aa->i; aj=aa->j; adiag=aa->diag; 74916ebf90aSShri Abhyankar bi=bb->i; bj=bb->j; garray = mat->garray; 75016ebf90aSShri Abhyankar av=aa->a; bv=bb->a; 7512205254eSKarl Rupp 75216ebf90aSShri Abhyankar rstart = A->rmap->rstart; 75316ebf90aSShri Abhyankar 754bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 755e0bace9bSHong Zhang nza = 0; /* num of upper triangular entries in mat->A, including diagonals */ 756e0bace9bSHong Zhang nzb = 0; /* num of upper triangular entries in mat->B */ 75716ebf90aSShri Abhyankar for (i=0; i<m; i++) { 758e0bace9bSHong Zhang nza += (ai[i+1] - adiag[i]); 75916ebf90aSShri Abhyankar countB = bi[i+1] - bi[i]; 76016ebf90aSShri Abhyankar bjj = bj + bi[i]; 761e0bace9bSHong Zhang for (j=0; j<countB; j++) { 762e0bace9bSHong Zhang if (garray[bjj[j]] > rstart) nzb++; 763e0bace9bSHong Zhang } 764e0bace9bSHong Zhang } 76516ebf90aSShri Abhyankar 766e0bace9bSHong Zhang nz = nza + nzb; /* total nz of upper triangular part of mat */ 76716ebf90aSShri Abhyankar *nnz = nz; 768185f6596SHong Zhang ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 769185f6596SHong Zhang col = row + nz; 770185f6596SHong Zhang val = (PetscScalar*)(col + nz); 771185f6596SHong Zhang 77216ebf90aSShri Abhyankar *r = row; *c = col; *v = val; 77316ebf90aSShri Abhyankar } else { 77416ebf90aSShri Abhyankar row = *r; col = *c; val = *v; 77516ebf90aSShri Abhyankar } 77616ebf90aSShri Abhyankar 77716ebf90aSShri Abhyankar jj = 0; irow = rstart; 77816ebf90aSShri Abhyankar for (i=0; i<m; i++) { 77916ebf90aSShri Abhyankar ajj = aj + adiag[i]; /* ptr to the beginning of the diagonal of this row */ 78016ebf90aSShri Abhyankar v1 = av + adiag[i]; 78116ebf90aSShri Abhyankar countA = ai[i+1] - adiag[i]; 78216ebf90aSShri Abhyankar countB = bi[i+1] - bi[i]; 78316ebf90aSShri Abhyankar bjj = bj + bi[i]; 78416ebf90aSShri Abhyankar v2 = bv + bi[i]; 78516ebf90aSShri Abhyankar 78616ebf90aSShri Abhyankar /* A-part */ 78716ebf90aSShri Abhyankar for (j=0; j<countA; j++) { 788bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 78916ebf90aSShri Abhyankar row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 79016ebf90aSShri Abhyankar } 79116ebf90aSShri Abhyankar val[jj++] = v1[j]; 79216ebf90aSShri Abhyankar } 79316ebf90aSShri Abhyankar 79416ebf90aSShri Abhyankar /* B-part */ 79516ebf90aSShri Abhyankar for (j=0; j < countB; j++) { 79616ebf90aSShri Abhyankar if (garray[bjj[j]] > rstart) { 797bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 79816ebf90aSShri Abhyankar row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 79916ebf90aSShri Abhyankar } 80016ebf90aSShri Abhyankar val[jj++] = v2[j]; 80116ebf90aSShri Abhyankar } 802397b6df1SKris Buschelman } 803397b6df1SKris Buschelman irow++; 804397b6df1SKris Buschelman } 805397b6df1SKris Buschelman PetscFunctionReturn(0); 806397b6df1SKris Buschelman } 807397b6df1SKris Buschelman 808397b6df1SKris Buschelman #undef __FUNCT__ 80920be8e61SHong Zhang #define __FUNCT__ "MatGetDiagonal_MUMPS" 81020be8e61SHong Zhang PetscErrorCode MatGetDiagonal_MUMPS(Mat A,Vec v) 81120be8e61SHong Zhang { 81220be8e61SHong Zhang PetscFunctionBegin; 81320be8e61SHong Zhang SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type: MUMPS factor"); 81420be8e61SHong Zhang PetscFunctionReturn(0); 81520be8e61SHong Zhang } 81620be8e61SHong Zhang 81720be8e61SHong Zhang #undef __FUNCT__ 8183924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS" 819dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A) 820dfbe8321SBarry Smith { 821a5e57a09SHong Zhang Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 822dfbe8321SBarry Smith PetscErrorCode ierr; 823b24902e0SBarry Smith 824397b6df1SKris Buschelman PetscFunctionBegin; 825a5e57a09SHong Zhang ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr); 826a5e57a09SHong Zhang ierr = VecScatterDestroy(&mumps->scat_rhs);CHKERRQ(ierr); 827a5e57a09SHong Zhang ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr); 828801fbe65SHong Zhang ierr = VecDestroy(&mumps->b_seq);CHKERRQ(ierr); 829a5e57a09SHong Zhang ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr); 830a5e57a09SHong Zhang ierr = PetscFree(mumps->id.perm_in);CHKERRQ(ierr); 831a5e57a09SHong Zhang ierr = PetscFree(mumps->irn);CHKERRQ(ierr); 832b34f08ffSHong Zhang ierr = PetscFree(mumps->info);CHKERRQ(ierr); 83359ac8732SStefano Zampini ierr = MatMumpsResetSchur_Private(mumps);CHKERRQ(ierr); 834a5e57a09SHong Zhang mumps->id.job = JOB_END; 835a5e57a09SHong Zhang PetscMUMPS_c(&mumps->id); 8366f3cc6f9SBarry Smith ierr = MPI_Comm_free(&mumps->comm_mumps);CHKERRQ(ierr); 837a5e57a09SHong Zhang if (mumps->Destroy) { 838a5e57a09SHong Zhang ierr = (mumps->Destroy)(A);CHKERRQ(ierr); 839bf0cc555SLisandro Dalcin } 840bf0cc555SLisandro Dalcin ierr = PetscFree(A->spptr);CHKERRQ(ierr); 841bf0cc555SLisandro Dalcin 84297969023SHong Zhang /* clear composed functions */ 843bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr); 8445a05ddb0SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorSetSchurIS_C",NULL);CHKERRQ(ierr); 8455a05ddb0SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorInvertSchurComplement_C",NULL);CHKERRQ(ierr); 8465a05ddb0SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorCreateSchurComplement_C",NULL);CHKERRQ(ierr); 8475a05ddb0SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSchurComplement_C",NULL);CHKERRQ(ierr); 8485a05ddb0SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorSolveSchurComplement_C",NULL);CHKERRQ(ierr); 8495a05ddb0SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorSolveSchurComplementTranspose_C",NULL);CHKERRQ(ierr); 850bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetIcntl_C",NULL);CHKERRQ(ierr); 851bc6112feSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetIcntl_C",NULL);CHKERRQ(ierr); 852bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetCntl_C",NULL);CHKERRQ(ierr); 853bc6112feSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetCntl_C",NULL);CHKERRQ(ierr); 854ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfo_C",NULL);CHKERRQ(ierr); 855ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfog_C",NULL);CHKERRQ(ierr); 856ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfo_C",NULL);CHKERRQ(ierr); 857ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfog_C",NULL);CHKERRQ(ierr); 858397b6df1SKris Buschelman PetscFunctionReturn(0); 859397b6df1SKris Buschelman } 860397b6df1SKris Buschelman 861397b6df1SKris Buschelman #undef __FUNCT__ 862f6c57405SHong Zhang #define __FUNCT__ "MatSolve_MUMPS" 863b24902e0SBarry Smith PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x) 864b24902e0SBarry Smith { 865a5e57a09SHong Zhang Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 866d54de34fSKris Buschelman PetscScalar *array; 86767877ebaSShri Abhyankar Vec b_seq; 868329ec9b3SHong Zhang IS is_iden,is_petsc; 869dfbe8321SBarry Smith PetscErrorCode ierr; 870329ec9b3SHong Zhang PetscInt i; 871cc86f929SStefano Zampini PetscBool second_solve = PETSC_FALSE; 872883f2eb9SBarry Smith static PetscBool cite1 = PETSC_FALSE,cite2 = PETSC_FALSE; 873397b6df1SKris Buschelman 874397b6df1SKris Buschelman PetscFunctionBegin; 875883f2eb9SBarry 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); 876883f2eb9SBarry 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); 8772aca8efcSHong Zhang 878d8cd2189SHong Zhang if (A->errortype) { 8792aca8efcSHong Zhang ierr = PetscInfo2(A,"MatSolve is called with singular matrix factor, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr); 8802aca8efcSHong Zhang ierr = VecSetInf(x);CHKERRQ(ierr); 8812aca8efcSHong Zhang PetscFunctionReturn(0); 8822aca8efcSHong Zhang } 8832aca8efcSHong Zhang 884a5e57a09SHong Zhang mumps->id.nrhs = 1; 885a5e57a09SHong Zhang b_seq = mumps->b_seq; 886a5e57a09SHong Zhang if (mumps->size > 1) { 887329ec9b3SHong Zhang /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */ 888a5e57a09SHong Zhang ierr = VecScatterBegin(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 889a5e57a09SHong Zhang ierr = VecScatterEnd(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 890a5e57a09SHong Zhang if (!mumps->myid) {ierr = VecGetArray(b_seq,&array);CHKERRQ(ierr);} 891397b6df1SKris Buschelman } else { /* size == 1 */ 892397b6df1SKris Buschelman ierr = VecCopy(b,x);CHKERRQ(ierr); 893397b6df1SKris Buschelman ierr = VecGetArray(x,&array);CHKERRQ(ierr); 894397b6df1SKris Buschelman } 895a5e57a09SHong Zhang if (!mumps->myid) { /* define rhs on the host */ 896a5e57a09SHong Zhang mumps->id.nrhs = 1; 897940cd9d6SSatish Balay mumps->id.rhs = (MumpsScalar*)array; 898397b6df1SKris Buschelman } 899397b6df1SKris Buschelman 900cc86f929SStefano Zampini /* 901cc86f929SStefano Zampini handle condensation step of Schur complement (if any) 902cc86f929SStefano Zampini We set by default ICNTL(26) == -1 when Schur indices have been provided by the user. 903cc86f929SStefano Zampini According to MUMPS (5.0.0) manual, any value should be harmful during the factorization phase 904cc86f929SStefano Zampini Unless the user provides a valid value for ICNTL(26), MatSolve and MatMatSolve routines solve the full system. 905cc86f929SStefano Zampini This requires an extra call to PetscMUMPS_c and the computation of the factors for S 906cc86f929SStefano Zampini */ 907cc86f929SStefano Zampini if (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2) { 908cc86f929SStefano Zampini second_solve = PETSC_TRUE; 909b8f61ee1SStefano Zampini ierr = MatMumpsHandleSchur_Private(mumps,PETSC_FALSE);CHKERRQ(ierr); 910cc86f929SStefano Zampini } 911397b6df1SKris Buschelman /* solve phase */ 912329ec9b3SHong Zhang /*-------------*/ 913a5e57a09SHong Zhang mumps->id.job = JOB_SOLVE; 914a5e57a09SHong Zhang PetscMUMPS_c(&mumps->id); 915a5e57a09SHong 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)); 916397b6df1SKris Buschelman 917b5fa320bSStefano Zampini /* handle expansion step of Schur complement (if any) */ 918cc86f929SStefano Zampini if (second_solve) { 919b8f61ee1SStefano Zampini ierr = MatMumpsHandleSchur_Private(mumps,PETSC_TRUE);CHKERRQ(ierr); 920cc86f929SStefano Zampini } 921b5fa320bSStefano Zampini 922a5e57a09SHong Zhang if (mumps->size > 1) { /* convert mumps distributed solution to petsc mpi x */ 923a5e57a09SHong Zhang if (mumps->scat_sol && mumps->ICNTL9_pre != mumps->id.ICNTL(9)) { 924a5e57a09SHong Zhang /* when id.ICNTL(9) changes, the contents of lsol_loc may change (not its size, lsol_loc), recreates scat_sol */ 925a5e57a09SHong Zhang ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr); 926397b6df1SKris Buschelman } 927a5e57a09SHong Zhang if (!mumps->scat_sol) { /* create scatter scat_sol */ 928a5e57a09SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,mumps->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */ 929a5e57a09SHong Zhang for (i=0; i<mumps->id.lsol_loc; i++) { 930a5e57a09SHong Zhang mumps->id.isol_loc[i] -= 1; /* change Fortran style to C style */ 931a5e57a09SHong Zhang } 932a5e57a09SHong Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,mumps->id.lsol_loc,mumps->id.isol_loc,PETSC_COPY_VALUES,&is_petsc);CHKERRQ(ierr); /* to */ 933a5e57a09SHong Zhang ierr = VecScatterCreate(mumps->x_seq,is_iden,x,is_petsc,&mumps->scat_sol);CHKERRQ(ierr); 9346bf464f9SBarry Smith ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 9356bf464f9SBarry Smith ierr = ISDestroy(&is_petsc);CHKERRQ(ierr); 9362205254eSKarl Rupp 937a5e57a09SHong Zhang mumps->ICNTL9_pre = mumps->id.ICNTL(9); /* save current value of id.ICNTL(9) */ 938397b6df1SKris Buschelman } 939a5e57a09SHong Zhang 940a5e57a09SHong Zhang ierr = VecScatterBegin(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 941a5e57a09SHong Zhang ierr = VecScatterEnd(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 942329ec9b3SHong Zhang } 943397b6df1SKris Buschelman PetscFunctionReturn(0); 944397b6df1SKris Buschelman } 945397b6df1SKris Buschelman 94651d5961aSHong Zhang #undef __FUNCT__ 94751d5961aSHong Zhang #define __FUNCT__ "MatSolveTranspose_MUMPS" 94851d5961aSHong Zhang PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x) 94951d5961aSHong Zhang { 950a5e57a09SHong Zhang Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 95151d5961aSHong Zhang PetscErrorCode ierr; 95251d5961aSHong Zhang 95351d5961aSHong Zhang PetscFunctionBegin; 954a5e57a09SHong Zhang mumps->id.ICNTL(9) = 0; 9550ad0caddSJed Brown ierr = MatSolve_MUMPS(A,b,x);CHKERRQ(ierr); 956a5e57a09SHong Zhang mumps->id.ICNTL(9) = 1; 95751d5961aSHong Zhang PetscFunctionReturn(0); 95851d5961aSHong Zhang } 95951d5961aSHong Zhang 960e0b74bf9SHong Zhang #undef __FUNCT__ 961e0b74bf9SHong Zhang #define __FUNCT__ "MatMatSolve_MUMPS" 962e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X) 963e0b74bf9SHong Zhang { 964bda8bf91SBarry Smith PetscErrorCode ierr; 965bda8bf91SBarry Smith PetscBool flg; 9664e34a73bSHong Zhang Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 967334c5f61SHong Zhang PetscInt i,nrhs,M; 9682cd7d884SHong Zhang PetscScalar *array,*bray; 969bda8bf91SBarry Smith 970e0b74bf9SHong Zhang PetscFunctionBegin; 9710298fd71SBarry Smith ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 972801fbe65SHong Zhang if (!flg) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 9730298fd71SBarry Smith ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 974801fbe65SHong Zhang if (!flg) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); 975801fbe65SHong 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"); 9764e34a73bSHong Zhang 9772cd7d884SHong Zhang ierr = MatGetSize(B,&M,&nrhs);CHKERRQ(ierr); 978334c5f61SHong Zhang mumps->id.nrhs = nrhs; 979334c5f61SHong Zhang mumps->id.lrhs = M; 9804e34a73bSHong Zhang 9812cd7d884SHong Zhang if (mumps->size == 1) { 9822cd7d884SHong Zhang /* copy B to X */ 9832cd7d884SHong Zhang ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr); 9842cd7d884SHong Zhang ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr); 9856444a565SStefano Zampini ierr = PetscMemcpy(array,bray,M*nrhs*sizeof(PetscScalar));CHKERRQ(ierr); 9862cd7d884SHong Zhang ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr); 987940cd9d6SSatish Balay mumps->id.rhs = (MumpsScalar*)array; 988b5fa320bSStefano Zampini /* handle condensation step of Schur complement (if any) */ 989b8f61ee1SStefano Zampini ierr = MatMumpsHandleSchur_Private(mumps,PETSC_FALSE);CHKERRQ(ierr); 990801fbe65SHong Zhang 9912cd7d884SHong Zhang /* solve phase */ 9922cd7d884SHong Zhang /*-------------*/ 9932cd7d884SHong Zhang mumps->id.job = JOB_SOLVE; 9942cd7d884SHong Zhang PetscMUMPS_c(&mumps->id); 9952cd7d884SHong 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)); 996b5fa320bSStefano Zampini 997b5fa320bSStefano Zampini /* handle expansion step of Schur complement (if any) */ 998b8f61ee1SStefano Zampini ierr = MatMumpsHandleSchur_Private(mumps,PETSC_TRUE);CHKERRQ(ierr); 9992cd7d884SHong Zhang ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr); 1000334c5f61SHong Zhang } else { /*--------- parallel case --------*/ 100171aed81dSHong Zhang PetscInt lsol_loc,nlsol_loc,*isol_loc,*idx,*iidx,*idxx,*isol_loc_save; 10021070efccSSatish Balay MumpsScalar *sol_loc,*sol_loc_save; 1003801fbe65SHong Zhang IS is_to,is_from; 1004334c5f61SHong Zhang PetscInt k,proc,j,m; 1005801fbe65SHong Zhang const PetscInt *rstart; 1006334c5f61SHong Zhang Vec v_mpi,b_seq,x_seq; 1007334c5f61SHong Zhang VecScatter scat_rhs,scat_sol; 1008801fbe65SHong Zhang 1009801fbe65SHong Zhang /* create x_seq to hold local solution */ 101071aed81dSHong Zhang isol_loc_save = mumps->id.isol_loc; /* save it for MatSovle() */ 101171aed81dSHong Zhang sol_loc_save = mumps->id.sol_loc; 1012801fbe65SHong Zhang 101371aed81dSHong Zhang lsol_loc = mumps->id.INFO(23); 101471aed81dSHong Zhang nlsol_loc = nrhs*lsol_loc; /* length of sol_loc */ 101571aed81dSHong Zhang ierr = PetscMalloc2(nlsol_loc,&sol_loc,nlsol_loc,&isol_loc);CHKERRQ(ierr); 1016940cd9d6SSatish Balay mumps->id.sol_loc = (MumpsScalar*)sol_loc; 1017801fbe65SHong Zhang mumps->id.isol_loc = isol_loc; 1018801fbe65SHong Zhang 10191070efccSSatish Balay ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,nlsol_loc,(PetscScalar*)sol_loc,&x_seq);CHKERRQ(ierr); 10202cd7d884SHong Zhang 102174f0fcc7SHong Zhang /* copy rhs matrix B into vector v_mpi */ 1022334c5f61SHong Zhang ierr = MatGetLocalSize(B,&m,NULL);CHKERRQ(ierr); 1023801fbe65SHong Zhang ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr); 102474f0fcc7SHong Zhang ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)B),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi);CHKERRQ(ierr); 1025801fbe65SHong Zhang ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr); 1026801fbe65SHong Zhang 1027334c5f61SHong Zhang /* scatter v_mpi to b_seq because MUMPS only supports centralized rhs */ 102874f0fcc7SHong Zhang /* idx: maps from k-th index of v_mpi to (i,j)-th global entry of B; 1029801fbe65SHong Zhang iidx: inverse of idx, will be used by scattering xx_seq -> X */ 1030801fbe65SHong Zhang ierr = PetscMalloc2(nrhs*M,&idx,nrhs*M,&iidx);CHKERRQ(ierr); 1031801fbe65SHong Zhang ierr = MatGetOwnershipRanges(B,&rstart);CHKERRQ(ierr); 1032801fbe65SHong Zhang k = 0; 1033801fbe65SHong Zhang for (proc=0; proc<mumps->size; proc++){ 1034801fbe65SHong Zhang for (j=0; j<nrhs; j++){ 1035801fbe65SHong Zhang for (i=rstart[proc]; i<rstart[proc+1]; i++){ 1036801fbe65SHong Zhang iidx[j*M + i] = k; 1037801fbe65SHong Zhang idx[k++] = j*M + i; 1038801fbe65SHong Zhang } 1039801fbe65SHong Zhang } 10402cd7d884SHong Zhang } 10412cd7d884SHong Zhang 1042801fbe65SHong Zhang if (!mumps->myid) { 1043334c5f61SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,nrhs*M,&b_seq);CHKERRQ(ierr); 1044801fbe65SHong Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,nrhs*M,idx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr); 1045801fbe65SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,nrhs*M,0,1,&is_from);CHKERRQ(ierr); 1046801fbe65SHong Zhang } else { 1047334c5f61SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,0,&b_seq);CHKERRQ(ierr); 1048801fbe65SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_to);CHKERRQ(ierr); 1049801fbe65SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_from);CHKERRQ(ierr); 1050801fbe65SHong Zhang } 1051334c5f61SHong Zhang ierr = VecScatterCreate(v_mpi,is_from,b_seq,is_to,&scat_rhs);CHKERRQ(ierr); 1052334c5f61SHong Zhang ierr = VecScatterBegin(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1053801fbe65SHong Zhang ierr = ISDestroy(&is_to);CHKERRQ(ierr); 1054801fbe65SHong Zhang ierr = ISDestroy(&is_from);CHKERRQ(ierr); 1055334c5f61SHong Zhang ierr = VecScatterEnd(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1056801fbe65SHong Zhang 1057801fbe65SHong Zhang if (!mumps->myid) { /* define rhs on the host */ 1058334c5f61SHong Zhang ierr = VecGetArray(b_seq,&bray);CHKERRQ(ierr); 1059940cd9d6SSatish Balay mumps->id.rhs = (MumpsScalar*)bray; 1060334c5f61SHong Zhang ierr = VecRestoreArray(b_seq,&bray);CHKERRQ(ierr); 1061801fbe65SHong Zhang } 1062801fbe65SHong Zhang 1063801fbe65SHong Zhang /* solve phase */ 1064801fbe65SHong Zhang /*-------------*/ 1065801fbe65SHong Zhang mumps->id.job = JOB_SOLVE; 1066801fbe65SHong Zhang PetscMUMPS_c(&mumps->id); 1067801fbe65SHong 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)); 1068801fbe65SHong Zhang 1069334c5f61SHong Zhang /* scatter mumps distributed solution to petsc vector v_mpi, which shares local arrays with solution matrix X */ 107074f0fcc7SHong Zhang ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr); 107174f0fcc7SHong Zhang ierr = VecPlaceArray(v_mpi,array);CHKERRQ(ierr); 1072801fbe65SHong Zhang 1073334c5f61SHong Zhang /* create scatter scat_sol */ 107471aed81dSHong Zhang ierr = PetscMalloc1(nlsol_loc,&idxx);CHKERRQ(ierr); 107571aed81dSHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,nlsol_loc,0,1,&is_from);CHKERRQ(ierr); 107671aed81dSHong Zhang for (i=0; i<lsol_loc; i++) { 1077334c5f61SHong Zhang isol_loc[i] -= 1; /* change Fortran style to C style */ 1078334c5f61SHong Zhang idxx[i] = iidx[isol_loc[i]]; 1079801fbe65SHong Zhang for (j=1; j<nrhs; j++){ 1080334c5f61SHong Zhang idxx[j*lsol_loc+i] = iidx[isol_loc[i]+j*M]; 1081801fbe65SHong Zhang } 1082801fbe65SHong Zhang } 108371aed81dSHong Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,nlsol_loc,idxx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr); 1084334c5f61SHong Zhang ierr = VecScatterCreate(x_seq,is_from,v_mpi,is_to,&scat_sol);CHKERRQ(ierr); 1085334c5f61SHong Zhang ierr = VecScatterBegin(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1086801fbe65SHong Zhang ierr = ISDestroy(&is_from);CHKERRQ(ierr); 1087801fbe65SHong Zhang ierr = ISDestroy(&is_to);CHKERRQ(ierr); 1088334c5f61SHong Zhang ierr = VecScatterEnd(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1089801fbe65SHong Zhang ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr); 109071aed81dSHong Zhang 109171aed81dSHong Zhang /* free spaces */ 109271aed81dSHong Zhang mumps->id.sol_loc = sol_loc_save; 109371aed81dSHong Zhang mumps->id.isol_loc = isol_loc_save; 109471aed81dSHong Zhang 109571aed81dSHong Zhang ierr = PetscFree2(sol_loc,isol_loc);CHKERRQ(ierr); 1096801fbe65SHong Zhang ierr = PetscFree2(idx,iidx);CHKERRQ(ierr); 1097801fbe65SHong Zhang ierr = PetscFree(idxx);CHKERRQ(ierr); 109871aed81dSHong Zhang ierr = VecDestroy(&x_seq);CHKERRQ(ierr); 109974f0fcc7SHong Zhang ierr = VecDestroy(&v_mpi);CHKERRQ(ierr); 1100334c5f61SHong Zhang ierr = VecDestroy(&b_seq);CHKERRQ(ierr); 1101334c5f61SHong Zhang ierr = VecScatterDestroy(&scat_rhs);CHKERRQ(ierr); 1102334c5f61SHong Zhang ierr = VecScatterDestroy(&scat_sol);CHKERRQ(ierr); 1103801fbe65SHong Zhang } 1104e0b74bf9SHong Zhang PetscFunctionReturn(0); 1105e0b74bf9SHong Zhang } 1106e0b74bf9SHong Zhang 1107ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX) 1108a58c3f20SHong Zhang /* 1109a58c3f20SHong Zhang input: 1110a58c3f20SHong Zhang F: numeric factor 1111a58c3f20SHong Zhang output: 1112a58c3f20SHong Zhang nneg: total number of negative pivots 1113a58c3f20SHong Zhang nzero: 0 1114a58c3f20SHong Zhang npos: (global dimension of F) - nneg 1115a58c3f20SHong Zhang */ 1116a58c3f20SHong Zhang 1117a58c3f20SHong Zhang #undef __FUNCT__ 1118a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS" 1119dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos) 1120a58c3f20SHong Zhang { 1121a5e57a09SHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 1122dfbe8321SBarry Smith PetscErrorCode ierr; 1123c1490034SHong Zhang PetscMPIInt size; 1124a58c3f20SHong Zhang 1125a58c3f20SHong Zhang PetscFunctionBegin; 1126ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)F),&size);CHKERRQ(ierr); 1127bcb30aebSHong 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 */ 1128a5e57a09SHong 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)); 1129ed85ac9fSHong Zhang 1130710ac8efSHong Zhang if (nneg) *nneg = mumps->id.INFOG(12); 1131ed85ac9fSHong Zhang if (nzero || npos) { 1132ed85ac9fSHong 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"); 1133710ac8efSHong Zhang if (nzero) *nzero = mumps->id.INFOG(28); 1134710ac8efSHong Zhang if (npos) *npos = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28)); 1135a58c3f20SHong Zhang } 1136a58c3f20SHong Zhang PetscFunctionReturn(0); 1137a58c3f20SHong Zhang } 1138ace3df97SHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */ 1139a58c3f20SHong Zhang 1140397b6df1SKris Buschelman #undef __FUNCT__ 1141f6c57405SHong Zhang #define __FUNCT__ "MatFactorNumeric_MUMPS" 11429a625307SHong Zhang PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info) 1143af281ebdSHong Zhang { 1144a5e57a09SHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)(F)->spptr; 11456849ba73SBarry Smith PetscErrorCode ierr; 1146e09efc27SHong Zhang Mat F_diag; 1147ace3abfcSBarry Smith PetscBool isMPIAIJ; 1148397b6df1SKris Buschelman 1149397b6df1SKris Buschelman PetscFunctionBegin; 11506baea169SHong Zhang if (mumps->id.INFOG(1) < 0) { 11512aca8efcSHong Zhang if (mumps->id.INFOG(1) == -6) { 11522aca8efcSHong Zhang ierr = PetscInfo2(A,"MatFactorNumeric is called with singular matrix structure, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr); 11536baea169SHong Zhang } 11546baea169SHong Zhang ierr = PetscInfo2(A,"MatFactorNumeric is called after analysis phase fails, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr); 11552aca8efcSHong Zhang PetscFunctionReturn(0); 11562aca8efcSHong Zhang } 11576baea169SHong Zhang 1158a5e57a09SHong Zhang ierr = (*mumps->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr); 1159397b6df1SKris Buschelman 1160397b6df1SKris Buschelman /* numerical factorization phase */ 1161329ec9b3SHong Zhang /*-------------------------------*/ 1162a5e57a09SHong Zhang mumps->id.job = JOB_FACTNUMERIC; 11634e34a73bSHong Zhang if (!mumps->id.ICNTL(18)) { /* A is centralized */ 1164a5e57a09SHong Zhang if (!mumps->myid) { 1165940cd9d6SSatish Balay mumps->id.a = (MumpsScalar*)mumps->val; 1166397b6df1SKris Buschelman } 1167397b6df1SKris Buschelman } else { 1168940cd9d6SSatish Balay mumps->id.a_loc = (MumpsScalar*)mumps->val; 1169397b6df1SKris Buschelman } 1170a5e57a09SHong Zhang PetscMUMPS_c(&mumps->id); 1171a5e57a09SHong Zhang if (mumps->id.INFOG(1) < 0) { 1172c0d63f2fSHong Zhang if (A->erroriffailure) { 1173c0d63f2fSHong Zhang SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2)); 1174c0d63f2fSHong Zhang } else { 1175c0d63f2fSHong Zhang if (mumps->id.INFOG(1) == -10) { /* numerically singular matrix */ 11762aca8efcSHong Zhang ierr = PetscInfo2(F,"matrix is numerically singular, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr); 1177680c5173SHong Zhang F->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 1178c0d63f2fSHong Zhang } else if (mumps->id.INFOG(1) == -13) { 1179c0d63f2fSHong Zhang ierr = PetscInfo2(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, cannot allocate required memory %d megabytes\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr); 1180680c5173SHong Zhang F->errortype = MAT_FACTOR_OUTMEMORY; 1181c0d63f2fSHong Zhang } else if (mumps->id.INFOG(1) == -8 || mumps->id.INFOG(1) == -9 || (-16 < mumps->id.INFOG(1) && mumps->id.INFOG(1) < -10) ) { 1182c0d63f2fSHong Zhang ierr = PetscInfo2(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d, problem with workarray \n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr); 1183680c5173SHong Zhang F->errortype = MAT_FACTOR_OUTMEMORY; 11842aca8efcSHong Zhang } else { 1185c0d63f2fSHong Zhang ierr = PetscInfo2(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr); 1186680c5173SHong Zhang F->errortype = MAT_FACTOR_OTHER; 11872aca8efcSHong Zhang } 11882aca8efcSHong Zhang } 1189397b6df1SKris Buschelman } 1190a5e57a09SHong 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)); 1191397b6df1SKris Buschelman 1192dcd589f8SShri Abhyankar (F)->assembled = PETSC_TRUE; 1193a5e57a09SHong Zhang mumps->matstruc = SAME_NONZERO_PATTERN; 1194b5fa320bSStefano Zampini mumps->schur_factored = PETSC_FALSE; 119559ac8732SStefano Zampini mumps->schur_inverted = PETSC_FALSE; 119667877ebaSShri Abhyankar 1197066565c5SStefano Zampini /* just to be sure that ICNTL(19) value returned by a call from MatMumpsGetIcntl is always consistent */ 1198066565c5SStefano Zampini if (!mumps->sym && mumps->id.ICNTL(19) && mumps->id.ICNTL(19) != 1) mumps->id.ICNTL(19) = 3; 1199066565c5SStefano Zampini 1200a5e57a09SHong Zhang if (mumps->size > 1) { 120167877ebaSShri Abhyankar PetscInt lsol_loc; 120267877ebaSShri Abhyankar PetscScalar *sol_loc; 12032205254eSKarl Rupp 1204c2093ab7SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr); 1205c2093ab7SHong Zhang if (isMPIAIJ) F_diag = ((Mat_MPIAIJ*)(F)->data)->A; 1206c2093ab7SHong Zhang else F_diag = ((Mat_MPISBAIJ*)(F)->data)->A; 1207c2093ab7SHong Zhang F_diag->assembled = PETSC_TRUE; 1208c2093ab7SHong Zhang 1209c2093ab7SHong Zhang /* distributed solution; Create x_seq=sol_loc for repeated use */ 1210c2093ab7SHong Zhang if (mumps->x_seq) { 1211c2093ab7SHong Zhang ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr); 1212c2093ab7SHong Zhang ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr); 1213c2093ab7SHong Zhang ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr); 1214c2093ab7SHong Zhang } 1215a5e57a09SHong Zhang lsol_loc = mumps->id.INFO(23); /* length of sol_loc */ 1216dcca6d9dSJed Brown ierr = PetscMalloc2(lsol_loc,&sol_loc,lsol_loc,&mumps->id.isol_loc);CHKERRQ(ierr); 1217a5e57a09SHong Zhang mumps->id.lsol_loc = lsol_loc; 1218940cd9d6SSatish Balay mumps->id.sol_loc = (MumpsScalar*)sol_loc; 1219a5e57a09SHong Zhang ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq);CHKERRQ(ierr); 122067877ebaSShri Abhyankar } 1221397b6df1SKris Buschelman PetscFunctionReturn(0); 1222397b6df1SKris Buschelman } 1223397b6df1SKris Buschelman 12249a2535b5SHong Zhang /* Sets MUMPS options from the options database */ 1225dcd589f8SShri Abhyankar #undef __FUNCT__ 12269a2535b5SHong Zhang #define __FUNCT__ "PetscSetMUMPSFromOptions" 12279a2535b5SHong Zhang PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A) 1228dcd589f8SShri Abhyankar { 12299a2535b5SHong Zhang Mat_MUMPS *mumps = (Mat_MUMPS*)F->spptr; 1230dcd589f8SShri Abhyankar PetscErrorCode ierr; 1231b34f08ffSHong Zhang PetscInt icntl,info[40],i,ninfo=40; 1232ace3abfcSBarry Smith PetscBool flg; 1233dcd589f8SShri Abhyankar 1234dcd589f8SShri Abhyankar PetscFunctionBegin; 1235ce94432eSBarry Smith ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr); 12369a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);CHKERRQ(ierr); 12379a2535b5SHong Zhang if (flg) mumps->id.ICNTL(1) = icntl; 12389a2535b5SHong 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); 12399a2535b5SHong Zhang if (flg) mumps->id.ICNTL(2) = icntl; 12409a2535b5SHong 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); 12419a2535b5SHong Zhang if (flg) mumps->id.ICNTL(3) = icntl; 1242dcd589f8SShri Abhyankar 12439a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",mumps->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr); 12449a2535b5SHong Zhang if (flg) mumps->id.ICNTL(4) = icntl; 12459a2535b5SHong Zhang if (mumps->id.ICNTL(4) || PetscLogPrintInfo) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */ 12469a2535b5SHong Zhang 1247d341cd04SHong 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); 12489a2535b5SHong Zhang if (flg) mumps->id.ICNTL(6) = icntl; 12499a2535b5SHong Zhang 1250d341cd04SHong 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); 1251dcd589f8SShri Abhyankar if (flg) { 12522205254eSKarl 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"); 12532205254eSKarl Rupp else mumps->id.ICNTL(7) = icntl; 1254dcd589f8SShri Abhyankar } 1255e0b74bf9SHong Zhang 12560298fd71SBarry 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); 1257d341cd04SHong 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() */ 12580298fd71SBarry 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); 1259d341cd04SHong 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); 1260d341cd04SHong 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); 1261d341cd04SHong 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); 1262d341cd04SHong 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); 1263d341cd04SHong 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); 126459ac8732SStefano Zampini if (mumps->id.ICNTL(19) <= 0 || mumps->id.ICNTL(19) > 3) { /* reset any schur data (if any) */ 126559ac8732SStefano Zampini ierr = MatMumpsResetSchur_Private(mumps);CHKERRQ(ierr); 126659ac8732SStefano Zampini } 12674e34a73bSHong 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 */ 1268d341cd04SHong 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 */ 12699a2535b5SHong Zhang 1270d341cd04SHong 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); 12710298fd71SBarry 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); 12720298fd71SBarry 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); 12739a2535b5SHong Zhang if (mumps->id.ICNTL(24)) { 12749a2535b5SHong Zhang mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */ 1275d7ebd59bSHong Zhang } 1276d7ebd59bSHong Zhang 1277d341cd04SHong 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); 1278d341cd04SHong 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); 12792cd7d884SHong 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); 12800298fd71SBarry 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); 1281d341cd04SHong 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); 12820298fd71SBarry 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); 1283d341cd04SHong 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); 12844e34a73bSHong 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 */ 12850298fd71SBarry Smith ierr = PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),NULL);CHKERRQ(ierr); 1286dcd589f8SShri Abhyankar 12870298fd71SBarry Smith ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),NULL);CHKERRQ(ierr); 12880298fd71SBarry 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); 12890298fd71SBarry Smith ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),NULL);CHKERRQ(ierr); 12900298fd71SBarry 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); 12910298fd71SBarry 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); 1292e5bb22a1SHong Zhang 12930298fd71SBarry Smith ierr = PetscOptionsString("-mat_mumps_ooc_tmpdir", "out of core directory", "None", mumps->id.ooc_tmpdir, mumps->id.ooc_tmpdir, 256, NULL); 1294b34f08ffSHong Zhang 129516d797efSHong Zhang ierr = PetscOptionsIntArray("-mat_mumps_view_info","request INFO local to each processor","",info,&ninfo,NULL);CHKERRQ(ierr); 1296b34f08ffSHong Zhang if (ninfo) { 1297b34f08ffSHong Zhang if (ninfo > 40) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"number of INFO %d must <= 40\n",ninfo); 1298b34f08ffSHong Zhang ierr = PetscMalloc1(ninfo,&mumps->info);CHKERRQ(ierr); 1299b34f08ffSHong Zhang mumps->ninfo = ninfo; 1300b34f08ffSHong Zhang for (i=0; i<ninfo; i++) { 13016c4ed002SBarry Smith if (info[i] < 0 || info[i]>40) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index of INFO %d must between 1 and 40\n",ninfo); 13026c4ed002SBarry Smith else { 1303b34f08ffSHong Zhang mumps->info[i] = info[i]; 1304b34f08ffSHong Zhang } 1305b34f08ffSHong Zhang } 1306b34f08ffSHong Zhang } 1307b34f08ffSHong Zhang 1308dcd589f8SShri Abhyankar PetscOptionsEnd(); 1309dcd589f8SShri Abhyankar PetscFunctionReturn(0); 1310dcd589f8SShri Abhyankar } 1311dcd589f8SShri Abhyankar 1312dcd589f8SShri Abhyankar #undef __FUNCT__ 1313dcd589f8SShri Abhyankar #define __FUNCT__ "PetscInitializeMUMPS" 1314f697e70eSHong Zhang PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps) 1315dcd589f8SShri Abhyankar { 1316dcd589f8SShri Abhyankar PetscErrorCode ierr; 1317dcd589f8SShri Abhyankar 1318dcd589f8SShri Abhyankar PetscFunctionBegin; 1319ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)A), &mumps->myid); 1320ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&mumps->size);CHKERRQ(ierr); 1321ce94432eSBarry Smith ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(mumps->comm_mumps));CHKERRQ(ierr); 13222205254eSKarl Rupp 1323f697e70eSHong Zhang mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps); 1324f697e70eSHong Zhang 1325f697e70eSHong Zhang mumps->id.job = JOB_INIT; 1326f697e70eSHong Zhang mumps->id.par = 1; /* host participates factorizaton and solve */ 1327f697e70eSHong Zhang mumps->id.sym = mumps->sym; 13282907cef9SHong Zhang PetscMUMPS_c(&mumps->id); 1329f697e70eSHong Zhang 13300298fd71SBarry Smith mumps->scat_rhs = NULL; 13310298fd71SBarry Smith mumps->scat_sol = NULL; 13329a2535b5SHong Zhang 133370544d5fSHong Zhang /* set PETSc-MUMPS default options - override MUMPS default */ 13349a2535b5SHong Zhang mumps->id.ICNTL(3) = 0; 13359a2535b5SHong Zhang mumps->id.ICNTL(4) = 0; 13369a2535b5SHong Zhang if (mumps->size == 1) { 13379a2535b5SHong Zhang mumps->id.ICNTL(18) = 0; /* centralized assembled matrix input */ 13389a2535b5SHong Zhang } else { 13399a2535b5SHong Zhang mumps->id.ICNTL(18) = 3; /* distributed assembled matrix input */ 13404e34a73bSHong Zhang mumps->id.ICNTL(20) = 0; /* rhs is in dense format */ 134170544d5fSHong Zhang mumps->id.ICNTL(21) = 1; /* distributed solution */ 13429a2535b5SHong Zhang } 13436444a565SStefano Zampini 13446444a565SStefano Zampini /* schur */ 13456444a565SStefano Zampini mumps->id.size_schur = 0; 13466444a565SStefano Zampini mumps->id.listvar_schur = NULL; 13476444a565SStefano Zampini mumps->id.schur = NULL; 1348b5fa320bSStefano Zampini mumps->sizeredrhs = 0; 1349b5fa320bSStefano Zampini mumps->schur_pivots = NULL; 1350b5fa320bSStefano Zampini mumps->schur_work = NULL; 135159ac8732SStefano Zampini mumps->schur_sol = NULL; 135259ac8732SStefano Zampini mumps->schur_sizesol = 0; 135359ac8732SStefano Zampini mumps->schur_factored = PETSC_FALSE; 135459ac8732SStefano Zampini mumps->schur_inverted = PETSC_FALSE; 1355dcd589f8SShri Abhyankar PetscFunctionReturn(0); 1356dcd589f8SShri Abhyankar } 1357dcd589f8SShri Abhyankar 13585cd7cf9dSHong Zhang #undef __FUNCT__ 13595cd7cf9dSHong Zhang #define __FUNCT__ "MatFactorSymbolic_MUMPS_ReportIfError" 13609a625307SHong Zhang PetscErrorCode MatFactorSymbolic_MUMPS_ReportIfError(Mat F,Mat A,const MatFactorInfo *info,Mat_MUMPS *mumps) 13615cd7cf9dSHong Zhang { 13625cd7cf9dSHong Zhang PetscErrorCode ierr; 13635cd7cf9dSHong Zhang 13645cd7cf9dSHong Zhang PetscFunctionBegin; 13655cd7cf9dSHong Zhang if (mumps->id.INFOG(1) < 0) { 13665cd7cf9dSHong Zhang if (A->erroriffailure) { 13675cd7cf9dSHong Zhang SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1)); 13685cd7cf9dSHong Zhang } else { 13695cd7cf9dSHong Zhang if (mumps->id.INFOG(1) == -6) { 13705cd7cf9dSHong Zhang ierr = PetscInfo2(F,"matrix is singular in structure, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr); 1371680c5173SHong Zhang F->errortype = MAT_FACTOR_STRUCT_ZEROPIVOT; 13725cd7cf9dSHong Zhang } else if (mumps->id.INFOG(1) == -5 || mumps->id.INFOG(1) == -7) { 13735cd7cf9dSHong Zhang ierr = PetscInfo2(F,"problem of workspace, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr); 1374680c5173SHong Zhang F->errortype = MAT_FACTOR_OUTMEMORY; 13755cd7cf9dSHong Zhang } else { 13765cd7cf9dSHong Zhang ierr = PetscInfo2(F,"Error reported by MUMPS in analysis phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr); 1377680c5173SHong Zhang F->errortype = MAT_FACTOR_OTHER; 13785cd7cf9dSHong Zhang } 13795cd7cf9dSHong Zhang } 13805cd7cf9dSHong Zhang } 13815cd7cf9dSHong Zhang PetscFunctionReturn(0); 13825cd7cf9dSHong Zhang } 13835cd7cf9dSHong Zhang 1384a5e57a09SHong 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 */ 1385397b6df1SKris Buschelman #undef __FUNCT__ 1386f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS" 13879a625307SHong Zhang PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 1388b24902e0SBarry Smith { 1389a5e57a09SHong Zhang Mat_MUMPS *mumps = (Mat_MUMPS*)F->spptr; 1390dcd589f8SShri Abhyankar PetscErrorCode ierr; 139167877ebaSShri Abhyankar Vec b; 139267877ebaSShri Abhyankar IS is_iden; 139367877ebaSShri Abhyankar const PetscInt M = A->rmap->N; 1394397b6df1SKris Buschelman 1395397b6df1SKris Buschelman PetscFunctionBegin; 1396a5e57a09SHong Zhang mumps->matstruc = DIFFERENT_NONZERO_PATTERN; 1397dcd589f8SShri Abhyankar 13989a2535b5SHong Zhang /* Set MUMPS options from the options database */ 13999a2535b5SHong Zhang ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr); 1400dcd589f8SShri Abhyankar 1401a5e57a09SHong Zhang ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr); 1402dcd589f8SShri Abhyankar 140367877ebaSShri Abhyankar /* analysis phase */ 140467877ebaSShri Abhyankar /*----------------*/ 1405a5e57a09SHong Zhang mumps->id.job = JOB_FACTSYMBOLIC; 1406a5e57a09SHong Zhang mumps->id.n = M; 1407a5e57a09SHong Zhang switch (mumps->id.ICNTL(18)) { 140867877ebaSShri Abhyankar case 0: /* centralized assembled matrix input */ 1409a5e57a09SHong Zhang if (!mumps->myid) { 1410a5e57a09SHong Zhang mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn; 1411a5e57a09SHong Zhang if (mumps->id.ICNTL(6)>1) { 1412940cd9d6SSatish Balay mumps->id.a = (MumpsScalar*)mumps->val; 141367877ebaSShri Abhyankar } 1414a5e57a09SHong Zhang if (mumps->id.ICNTL(7) == 1) { /* use user-provide matrix ordering - assuming r = c ordering */ 14155248a706SHong Zhang /* 14165248a706SHong Zhang PetscBool flag; 14175248a706SHong Zhang ierr = ISEqual(r,c,&flag);CHKERRQ(ierr); 14185248a706SHong Zhang if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"row_perm != col_perm"); 14195248a706SHong Zhang ierr = ISView(r,PETSC_VIEWER_STDOUT_SELF); 14205248a706SHong Zhang */ 1421a5e57a09SHong Zhang if (!mumps->myid) { 1422e0b74bf9SHong Zhang const PetscInt *idx; 1423e0b74bf9SHong Zhang PetscInt i,*perm_in; 14242205254eSKarl Rupp 1425785e854fSJed Brown ierr = PetscMalloc1(M,&perm_in);CHKERRQ(ierr); 1426e0b74bf9SHong Zhang ierr = ISGetIndices(r,&idx);CHKERRQ(ierr); 14272205254eSKarl Rupp 1428a5e57a09SHong Zhang mumps->id.perm_in = perm_in; 1429e0b74bf9SHong Zhang for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */ 1430e0b74bf9SHong Zhang ierr = ISRestoreIndices(r,&idx);CHKERRQ(ierr); 1431e0b74bf9SHong Zhang } 1432e0b74bf9SHong Zhang } 143367877ebaSShri Abhyankar } 143467877ebaSShri Abhyankar break; 143567877ebaSShri Abhyankar case 3: /* distributed assembled matrix input (size>1) */ 1436a5e57a09SHong Zhang mumps->id.nz_loc = mumps->nz; 1437a5e57a09SHong Zhang mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn; 1438a5e57a09SHong Zhang if (mumps->id.ICNTL(6)>1) { 1439940cd9d6SSatish Balay mumps->id.a_loc = (MumpsScalar*)mumps->val; 144067877ebaSShri Abhyankar } 144167877ebaSShri Abhyankar /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 1442a5e57a09SHong Zhang if (!mumps->myid) { 14432cd7d884SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,A->rmap->N,&mumps->b_seq);CHKERRQ(ierr); 14442cd7d884SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->N,0,1,&is_iden);CHKERRQ(ierr); 144567877ebaSShri Abhyankar } else { 1446a5e57a09SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr); 144767877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 144867877ebaSShri Abhyankar } 14492a7a6963SBarry Smith ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr); 1450a5e57a09SHong Zhang ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr); 14516bf464f9SBarry Smith ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 14526bf464f9SBarry Smith ierr = VecDestroy(&b);CHKERRQ(ierr); 145367877ebaSShri Abhyankar break; 145467877ebaSShri Abhyankar } 1455a5e57a09SHong Zhang PetscMUMPS_c(&mumps->id); 14565cd7cf9dSHong Zhang ierr = MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);CHKERRQ(ierr); 145767877ebaSShri Abhyankar 1458719d5645SBarry Smith F->ops->lufactornumeric = MatFactorNumeric_MUMPS; 1459dcd589f8SShri Abhyankar F->ops->solve = MatSolve_MUMPS; 146051d5961aSHong Zhang F->ops->solvetranspose = MatSolveTranspose_MUMPS; 14614e34a73bSHong Zhang F->ops->matsolve = MatMatSolve_MUMPS; 1462b24902e0SBarry Smith PetscFunctionReturn(0); 1463b24902e0SBarry Smith } 1464b24902e0SBarry Smith 1465450b117fSShri Abhyankar /* Note the Petsc r and c permutations are ignored */ 1466450b117fSShri Abhyankar #undef __FUNCT__ 1467450b117fSShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_BAIJMUMPS" 14689a625307SHong Zhang PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 1469450b117fSShri Abhyankar { 1470a5e57a09SHong Zhang Mat_MUMPS *mumps = (Mat_MUMPS*)F->spptr; 1471dcd589f8SShri Abhyankar PetscErrorCode ierr; 147267877ebaSShri Abhyankar Vec b; 147367877ebaSShri Abhyankar IS is_iden; 147467877ebaSShri Abhyankar const PetscInt M = A->rmap->N; 1475450b117fSShri Abhyankar 1476450b117fSShri Abhyankar PetscFunctionBegin; 1477a5e57a09SHong Zhang mumps->matstruc = DIFFERENT_NONZERO_PATTERN; 1478dcd589f8SShri Abhyankar 14799a2535b5SHong Zhang /* Set MUMPS options from the options database */ 14809a2535b5SHong Zhang ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr); 1481dcd589f8SShri Abhyankar 1482a5e57a09SHong Zhang ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr); 148367877ebaSShri Abhyankar 148467877ebaSShri Abhyankar /* analysis phase */ 148567877ebaSShri Abhyankar /*----------------*/ 1486a5e57a09SHong Zhang mumps->id.job = JOB_FACTSYMBOLIC; 1487a5e57a09SHong Zhang mumps->id.n = M; 1488a5e57a09SHong Zhang switch (mumps->id.ICNTL(18)) { 148967877ebaSShri Abhyankar case 0: /* centralized assembled matrix input */ 1490a5e57a09SHong Zhang if (!mumps->myid) { 1491a5e57a09SHong Zhang mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn; 1492a5e57a09SHong Zhang if (mumps->id.ICNTL(6)>1) { 1493940cd9d6SSatish Balay mumps->id.a = (MumpsScalar*)mumps->val; 149467877ebaSShri Abhyankar } 149567877ebaSShri Abhyankar } 149667877ebaSShri Abhyankar break; 149767877ebaSShri Abhyankar case 3: /* distributed assembled matrix input (size>1) */ 1498a5e57a09SHong Zhang mumps->id.nz_loc = mumps->nz; 1499a5e57a09SHong Zhang mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn; 1500a5e57a09SHong Zhang if (mumps->id.ICNTL(6)>1) { 1501940cd9d6SSatish Balay mumps->id.a_loc = (MumpsScalar*)mumps->val; 150267877ebaSShri Abhyankar } 150367877ebaSShri Abhyankar /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 1504a5e57a09SHong Zhang if (!mumps->myid) { 1505a5e57a09SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr); 150667877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 150767877ebaSShri Abhyankar } else { 1508a5e57a09SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr); 150967877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 151067877ebaSShri Abhyankar } 15112a7a6963SBarry Smith ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr); 1512a5e57a09SHong Zhang ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr); 15136bf464f9SBarry Smith ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 15146bf464f9SBarry Smith ierr = VecDestroy(&b);CHKERRQ(ierr); 151567877ebaSShri Abhyankar break; 151667877ebaSShri Abhyankar } 1517a5e57a09SHong Zhang PetscMUMPS_c(&mumps->id); 15185cd7cf9dSHong Zhang ierr = MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);CHKERRQ(ierr); 151967877ebaSShri Abhyankar 1520450b117fSShri Abhyankar F->ops->lufactornumeric = MatFactorNumeric_MUMPS; 1521dcd589f8SShri Abhyankar F->ops->solve = MatSolve_MUMPS; 152251d5961aSHong Zhang F->ops->solvetranspose = MatSolveTranspose_MUMPS; 1523450b117fSShri Abhyankar PetscFunctionReturn(0); 1524450b117fSShri Abhyankar } 1525b24902e0SBarry Smith 1526141f4205SHong Zhang /* Note the Petsc r permutation and factor info are ignored */ 1527397b6df1SKris Buschelman #undef __FUNCT__ 152867877ebaSShri Abhyankar #define __FUNCT__ "MatCholeskyFactorSymbolic_MUMPS" 15299a625307SHong Zhang PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info) 1530b24902e0SBarry Smith { 1531a5e57a09SHong Zhang Mat_MUMPS *mumps = (Mat_MUMPS*)F->spptr; 1532dcd589f8SShri Abhyankar PetscErrorCode ierr; 153367877ebaSShri Abhyankar Vec b; 153467877ebaSShri Abhyankar IS is_iden; 153567877ebaSShri Abhyankar const PetscInt M = A->rmap->N; 1536397b6df1SKris Buschelman 1537397b6df1SKris Buschelman PetscFunctionBegin; 1538a5e57a09SHong Zhang mumps->matstruc = DIFFERENT_NONZERO_PATTERN; 1539dcd589f8SShri Abhyankar 15409a2535b5SHong Zhang /* Set MUMPS options from the options database */ 15419a2535b5SHong Zhang ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr); 1542dcd589f8SShri Abhyankar 1543a5e57a09SHong Zhang ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr); 1544dcd589f8SShri Abhyankar 154567877ebaSShri Abhyankar /* analysis phase */ 154667877ebaSShri Abhyankar /*----------------*/ 1547a5e57a09SHong Zhang mumps->id.job = JOB_FACTSYMBOLIC; 1548a5e57a09SHong Zhang mumps->id.n = M; 1549a5e57a09SHong Zhang switch (mumps->id.ICNTL(18)) { 155067877ebaSShri Abhyankar case 0: /* centralized assembled matrix input */ 1551a5e57a09SHong Zhang if (!mumps->myid) { 1552a5e57a09SHong Zhang mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn; 1553a5e57a09SHong Zhang if (mumps->id.ICNTL(6)>1) { 1554940cd9d6SSatish Balay mumps->id.a = (MumpsScalar*)mumps->val; 155567877ebaSShri Abhyankar } 155667877ebaSShri Abhyankar } 155767877ebaSShri Abhyankar break; 155867877ebaSShri Abhyankar case 3: /* distributed assembled matrix input (size>1) */ 1559a5e57a09SHong Zhang mumps->id.nz_loc = mumps->nz; 1560a5e57a09SHong Zhang mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn; 1561a5e57a09SHong Zhang if (mumps->id.ICNTL(6)>1) { 1562940cd9d6SSatish Balay mumps->id.a_loc = (MumpsScalar*)mumps->val; 156367877ebaSShri Abhyankar } 156467877ebaSShri Abhyankar /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 1565a5e57a09SHong Zhang if (!mumps->myid) { 1566a5e57a09SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr); 156767877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 156867877ebaSShri Abhyankar } else { 1569a5e57a09SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr); 157067877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 157167877ebaSShri Abhyankar } 15722a7a6963SBarry Smith ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr); 1573a5e57a09SHong Zhang ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr); 15746bf464f9SBarry Smith ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 15756bf464f9SBarry Smith ierr = VecDestroy(&b);CHKERRQ(ierr); 157667877ebaSShri Abhyankar break; 157767877ebaSShri Abhyankar } 1578a5e57a09SHong Zhang PetscMUMPS_c(&mumps->id); 15795cd7cf9dSHong Zhang ierr = MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);CHKERRQ(ierr); 15805cd7cf9dSHong Zhang 158167877ebaSShri Abhyankar 15822792810eSHong Zhang F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS; 1583dcd589f8SShri Abhyankar F->ops->solve = MatSolve_MUMPS; 158451d5961aSHong Zhang F->ops->solvetranspose = MatSolve_MUMPS; 15854e34a73bSHong Zhang F->ops->matsolve = MatMatSolve_MUMPS; 15864e34a73bSHong Zhang #if defined(PETSC_USE_COMPLEX) 15870298fd71SBarry Smith F->ops->getinertia = NULL; 15884e34a73bSHong Zhang #else 15894e34a73bSHong Zhang F->ops->getinertia = MatGetInertia_SBAIJMUMPS; 1590db4efbfdSBarry Smith #endif 1591b24902e0SBarry Smith PetscFunctionReturn(0); 1592b24902e0SBarry Smith } 1593b24902e0SBarry Smith 1594397b6df1SKris Buschelman #undef __FUNCT__ 159564e6c443SBarry Smith #define __FUNCT__ "MatView_MUMPS" 159664e6c443SBarry Smith PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer) 159774ed9c26SBarry Smith { 1598f6c57405SHong Zhang PetscErrorCode ierr; 159964e6c443SBarry Smith PetscBool iascii; 160064e6c443SBarry Smith PetscViewerFormat format; 1601a5e57a09SHong Zhang Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 1602f6c57405SHong Zhang 1603f6c57405SHong Zhang PetscFunctionBegin; 160464e6c443SBarry Smith /* check if matrix is mumps type */ 160564e6c443SBarry Smith if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0); 160664e6c443SBarry Smith 1607251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 160864e6c443SBarry Smith if (iascii) { 160964e6c443SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 161064e6c443SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO) { 161164e6c443SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr); 1612a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," SYM (matrix type): %d \n",mumps->id.sym);CHKERRQ(ierr); 1613a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," PAR (host participation): %d \n",mumps->id.par);CHKERRQ(ierr); 1614a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(1) (output for error): %d \n",mumps->id.ICNTL(1));CHKERRQ(ierr); 1615a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(2) (output of diagnostic msg): %d \n",mumps->id.ICNTL(2));CHKERRQ(ierr); 1616a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(3) (output for global info): %d \n",mumps->id.ICNTL(3));CHKERRQ(ierr); 1617a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(4) (level of printing): %d \n",mumps->id.ICNTL(4));CHKERRQ(ierr); 1618a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(5) (input mat struct): %d \n",mumps->id.ICNTL(5));CHKERRQ(ierr); 1619a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(6) (matrix prescaling): %d \n",mumps->id.ICNTL(6));CHKERRQ(ierr); 1620a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(7) (sequentia matrix ordering):%d \n",mumps->id.ICNTL(7));CHKERRQ(ierr); 1621a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(8) (scalling strategy): %d \n",mumps->id.ICNTL(8));CHKERRQ(ierr); 1622a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(10) (max num of refinements): %d \n",mumps->id.ICNTL(10));CHKERRQ(ierr); 1623a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(11) (error analysis): %d \n",mumps->id.ICNTL(11));CHKERRQ(ierr); 1624a5e57a09SHong Zhang if (mumps->id.ICNTL(11)>0) { 1625a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(4) (inf norm of input mat): %g\n",mumps->id.RINFOG(4));CHKERRQ(ierr); 1626a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(5) (inf norm of solution): %g\n",mumps->id.RINFOG(5));CHKERRQ(ierr); 1627a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(6) (inf norm of residual): %g\n",mumps->id.RINFOG(6));CHKERRQ(ierr); 1628a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8));CHKERRQ(ierr); 1629a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(9) (error estimate): %g \n",mumps->id.RINFOG(9));CHKERRQ(ierr); 1630a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11));CHKERRQ(ierr); 1631f6c57405SHong Zhang } 1632a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(12) (efficiency control): %d \n",mumps->id.ICNTL(12));CHKERRQ(ierr); 1633a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(13) (efficiency control): %d \n",mumps->id.ICNTL(13));CHKERRQ(ierr); 1634a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(14) (percentage of estimated workspace increase): %d \n",mumps->id.ICNTL(14));CHKERRQ(ierr); 1635f6c57405SHong Zhang /* ICNTL(15-17) not used */ 1636a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(18) (input mat struct): %d \n",mumps->id.ICNTL(18));CHKERRQ(ierr); 1637a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(19) (Shur complement info): %d \n",mumps->id.ICNTL(19));CHKERRQ(ierr); 1638a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(20) (rhs sparse pattern): %d \n",mumps->id.ICNTL(20));CHKERRQ(ierr); 1639ca3dc52bSPierre Jolivet ierr = PetscViewerASCIIPrintf(viewer," ICNTL(21) (solution struct): %d \n",mumps->id.ICNTL(21));CHKERRQ(ierr); 1640a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(22) (in-core/out-of-core facility): %d \n",mumps->id.ICNTL(22));CHKERRQ(ierr); 1641a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(23) (max size of memory can be allocated locally):%d \n",mumps->id.ICNTL(23));CHKERRQ(ierr); 1642c0165424SHong Zhang 1643a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(24) (detection of null pivot rows): %d \n",mumps->id.ICNTL(24));CHKERRQ(ierr); 1644a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(25) (computation of a null space basis): %d \n",mumps->id.ICNTL(25));CHKERRQ(ierr); 1645a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(26) (Schur options for rhs or solution): %d \n",mumps->id.ICNTL(26));CHKERRQ(ierr); 1646a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(27) (experimental parameter): %d \n",mumps->id.ICNTL(27));CHKERRQ(ierr); 1647a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(28) (use parallel or sequential ordering): %d \n",mumps->id.ICNTL(28));CHKERRQ(ierr); 1648a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(29) (parallel ordering): %d \n",mumps->id.ICNTL(29));CHKERRQ(ierr); 164942179a6aSHong Zhang 1650a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(30) (user-specified set of entries in inv(A)): %d \n",mumps->id.ICNTL(30));CHKERRQ(ierr); 1651a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(31) (factors is discarded in the solve phase): %d \n",mumps->id.ICNTL(31));CHKERRQ(ierr); 1652a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(33) (compute determinant): %d \n",mumps->id.ICNTL(33));CHKERRQ(ierr); 1653f6c57405SHong Zhang 1654a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(1) (relative pivoting threshold): %g \n",mumps->id.CNTL(1));CHKERRQ(ierr); 1655a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2));CHKERRQ(ierr); 1656ca3dc52bSPierre Jolivet ierr = PetscViewerASCIIPrintf(viewer," CNTL(3) (absolute pivoting threshold): %g \n",mumps->id.CNTL(3));CHKERRQ(ierr); 1657ca3dc52bSPierre Jolivet ierr = PetscViewerASCIIPrintf(viewer," CNTL(4) (value of static pivoting): %g \n",mumps->id.CNTL(4));CHKERRQ(ierr); 1658a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(5) (fixation for null pivots): %g \n",mumps->id.CNTL(5));CHKERRQ(ierr); 1659f6c57405SHong Zhang 1660f6c57405SHong Zhang /* infomation local to each processor */ 166134ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr); 16621575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 1663a5e57a09SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",mumps->myid,mumps->id.RINFO(1));CHKERRQ(ierr); 166434ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 166534ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr); 1666a5e57a09SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",mumps->myid,mumps->id.RINFO(2));CHKERRQ(ierr); 166734ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 166834ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr); 1669a5e57a09SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",mumps->myid,mumps->id.RINFO(3));CHKERRQ(ierr); 167034ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 1671f6c57405SHong Zhang 167234ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr); 1673a5e57a09SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(15));CHKERRQ(ierr); 167434ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 1675f6c57405SHong Zhang 167634ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr); 1677a5e57a09SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(16));CHKERRQ(ierr); 167834ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 1679f6c57405SHong Zhang 168034ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr); 1681a5e57a09SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(23));CHKERRQ(ierr); 168234ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 1683b34f08ffSHong Zhang 1684b34f08ffSHong Zhang if (mumps->ninfo && mumps->ninfo <= 40){ 1685b34f08ffSHong Zhang PetscInt i; 1686b34f08ffSHong Zhang for (i=0; i<mumps->ninfo; i++){ 1687b34f08ffSHong Zhang ierr = PetscViewerASCIIPrintf(viewer, " INFO(%d): \n",mumps->info[i]);CHKERRQ(ierr); 1688b34f08ffSHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(mumps->info[i]));CHKERRQ(ierr); 1689b34f08ffSHong Zhang ierr = PetscViewerFlush(viewer); 1690b34f08ffSHong Zhang } 1691b34f08ffSHong Zhang } 1692b34f08ffSHong Zhang 1693b34f08ffSHong Zhang 16941575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 1695f6c57405SHong Zhang 1696a5e57a09SHong Zhang if (!mumps->myid) { /* information from the host */ 1697a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",mumps->id.RINFOG(1));CHKERRQ(ierr); 1698a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",mumps->id.RINFOG(2));CHKERRQ(ierr); 1699a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",mumps->id.RINFOG(3));CHKERRQ(ierr); 1700a5e57a09SHong 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); 1701f6c57405SHong Zhang 1702a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(3));CHKERRQ(ierr); 1703a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(4));CHKERRQ(ierr); 1704a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(5) (estimated maximum front size in the complete tree): %d \n",mumps->id.INFOG(5));CHKERRQ(ierr); 1705a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(6) (number of nodes in the complete tree): %d \n",mumps->id.INFOG(6));CHKERRQ(ierr); 1706a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(7) (ordering option effectively use after analysis): %d \n",mumps->id.INFOG(7));CHKERRQ(ierr); 1707a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",mumps->id.INFOG(8));CHKERRQ(ierr); 1708a5e57a09SHong 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); 1709a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(10) (total integer space store the matrix factors after factorization): %d \n",mumps->id.INFOG(10));CHKERRQ(ierr); 1710a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(11) (order of largest frontal matrix after factorization): %d \n",mumps->id.INFOG(11));CHKERRQ(ierr); 1711a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(12) (number of off-diagonal pivots): %d \n",mumps->id.INFOG(12));CHKERRQ(ierr); 1712a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(13) (number of delayed pivots after factorization): %d \n",mumps->id.INFOG(13));CHKERRQ(ierr); 1713a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(14) (number of memory compress after factorization): %d \n",mumps->id.INFOG(14));CHKERRQ(ierr); 1714a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(15) (number of steps of iterative refinement after solution): %d \n",mumps->id.INFOG(15));CHKERRQ(ierr); 1715a5e57a09SHong 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); 1716a5e57a09SHong 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); 1717a5e57a09SHong 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); 1718a5e57a09SHong 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); 1719a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(20) (estimated number of entries in the factors): %d \n",mumps->id.INFOG(20));CHKERRQ(ierr); 1720a5e57a09SHong 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); 1721a5e57a09SHong 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); 1722a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",mumps->id.INFOG(23));CHKERRQ(ierr); 1723a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",mumps->id.INFOG(24));CHKERRQ(ierr); 1724a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",mumps->id.INFOG(25));CHKERRQ(ierr); 172540d435e3SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(28) (after factorization: number of null pivots encountered): %d\n",mumps->id.INFOG(28));CHKERRQ(ierr); 172640d435e3SHong 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); 172740d435e3SHong 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); 172840d435e3SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(32) (after analysis: type of analysis done): %d\n",mumps->id.INFOG(32));CHKERRQ(ierr); 172940d435e3SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(33) (value used for ICNTL(8)): %d\n",mumps->id.INFOG(33));CHKERRQ(ierr); 173040d435e3SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(34) (exponent of the determinant if determinant is requested): %d\n",mumps->id.INFOG(34));CHKERRQ(ierr); 1731f6c57405SHong Zhang } 1732f6c57405SHong Zhang } 1733cb828f0fSHong Zhang } 1734f6c57405SHong Zhang PetscFunctionReturn(0); 1735f6c57405SHong Zhang } 1736f6c57405SHong Zhang 173735bd34faSBarry Smith #undef __FUNCT__ 173835bd34faSBarry Smith #define __FUNCT__ "MatGetInfo_MUMPS" 173935bd34faSBarry Smith PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info) 174035bd34faSBarry Smith { 1741cb828f0fSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)A->spptr; 174235bd34faSBarry Smith 174335bd34faSBarry Smith PetscFunctionBegin; 174435bd34faSBarry Smith info->block_size = 1.0; 1745cb828f0fSHong Zhang info->nz_allocated = mumps->id.INFOG(20); 1746cb828f0fSHong Zhang info->nz_used = mumps->id.INFOG(20); 174735bd34faSBarry Smith info->nz_unneeded = 0.0; 174835bd34faSBarry Smith info->assemblies = 0.0; 174935bd34faSBarry Smith info->mallocs = 0.0; 175035bd34faSBarry Smith info->memory = 0.0; 175135bd34faSBarry Smith info->fill_ratio_given = 0; 175235bd34faSBarry Smith info->fill_ratio_needed = 0; 175335bd34faSBarry Smith info->factor_mallocs = 0; 175435bd34faSBarry Smith PetscFunctionReturn(0); 175535bd34faSBarry Smith } 175635bd34faSBarry Smith 17575ccb76cbSHong Zhang /* -------------------------------------------------------------------------------------------*/ 17585ccb76cbSHong Zhang #undef __FUNCT__ 17598e7ba810SStefano Zampini #define __FUNCT__ "MatFactorSetSchurIS_MUMPS" 17608e7ba810SStefano Zampini PetscErrorCode MatFactorSetSchurIS_MUMPS(Mat F, IS is) 17616444a565SStefano Zampini { 17626444a565SStefano Zampini Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 17638e7ba810SStefano Zampini const PetscInt *idxs; 17648e7ba810SStefano Zampini PetscInt size,i; 17656444a565SStefano Zampini PetscErrorCode ierr; 17666444a565SStefano Zampini 17676444a565SStefano Zampini PetscFunctionBegin; 176859ac8732SStefano Zampini if (mumps->size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MUMPS parallel Schur complements not yet supported from PETSc\n"); 17698e7ba810SStefano Zampini ierr = ISGetLocalSize(is,&size);CHKERRQ(ierr); 17706444a565SStefano Zampini if (mumps->id.size_schur != size) { 17716444a565SStefano Zampini ierr = PetscFree2(mumps->id.listvar_schur,mumps->id.schur);CHKERRQ(ierr); 17726444a565SStefano Zampini mumps->id.size_schur = size; 17736444a565SStefano Zampini mumps->id.schur_lld = size; 17746444a565SStefano Zampini ierr = PetscMalloc2(size,&mumps->id.listvar_schur,size*size,&mumps->id.schur);CHKERRQ(ierr); 17756444a565SStefano Zampini } 17768e7ba810SStefano Zampini ierr = ISGetIndices(is,&idxs);CHKERRQ(ierr); 17776444a565SStefano Zampini ierr = PetscMemcpy(mumps->id.listvar_schur,idxs,size*sizeof(PetscInt));CHKERRQ(ierr); 17788e7ba810SStefano Zampini /* MUMPS expects Fortran style indices */ 17798e7ba810SStefano Zampini for (i=0;i<size;i++) mumps->id.listvar_schur[i]++; 17808e7ba810SStefano Zampini ierr = ISRestoreIndices(is,&idxs);CHKERRQ(ierr); 178178a7176cSStefano Zampini if (size) { /* turn on Schur switch if we the set of indices is not empty */ 17826444a565SStefano Zampini if (F->factortype == MAT_FACTOR_LU) { 178359ac8732SStefano Zampini mumps->id.ICNTL(19) = 3; /* MUMPS returns full matrix */ 17846444a565SStefano Zampini } else { 178559ac8732SStefano Zampini mumps->id.ICNTL(19) = 2; /* MUMPS returns lower triangular part */ 17866444a565SStefano Zampini } 178759ac8732SStefano Zampini /* set a special value of ICNTL (not handled my MUMPS) to be used in the solve phase by PETSc */ 1788b5fa320bSStefano Zampini mumps->id.ICNTL(26) = -1; 178978a7176cSStefano Zampini } 17906444a565SStefano Zampini PetscFunctionReturn(0); 17916444a565SStefano Zampini } 17926444a565SStefano Zampini 17936444a565SStefano Zampini /* -------------------------------------------------------------------------------------------*/ 17946444a565SStefano Zampini #undef __FUNCT__ 17955a05ddb0SStefano Zampini #define __FUNCT__ "MatFactorCreateSchurComplement_MUMPS" 17965a05ddb0SStefano Zampini PetscErrorCode MatFactorCreateSchurComplement_MUMPS(Mat F,Mat* S) 17976444a565SStefano Zampini { 17986444a565SStefano Zampini Mat St; 17996444a565SStefano Zampini Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 18006444a565SStefano Zampini PetscScalar *array; 18016444a565SStefano Zampini #if defined(PETSC_USE_COMPLEX) 18028ac429a0SStefano Zampini PetscScalar im = PetscSqrtScalar((PetscScalar)-1.0); 18036444a565SStefano Zampini #endif 18046444a565SStefano Zampini PetscErrorCode ierr; 18056444a565SStefano Zampini 18066444a565SStefano Zampini PetscFunctionBegin; 18075a05ddb0SStefano 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"); 18085a05ddb0SStefano Zampini else if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before"); 18096f3cc6f9SBarry Smith 18106444a565SStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)F),&St);CHKERRQ(ierr); 18116444a565SStefano Zampini ierr = MatSetSizes(St,PETSC_DECIDE,PETSC_DECIDE,mumps->id.size_schur,mumps->id.size_schur);CHKERRQ(ierr); 18126444a565SStefano Zampini ierr = MatSetType(St,MATDENSE);CHKERRQ(ierr); 18136444a565SStefano Zampini ierr = MatSetUp(St);CHKERRQ(ierr); 18146444a565SStefano Zampini ierr = MatDenseGetArray(St,&array);CHKERRQ(ierr); 181559ac8732SStefano Zampini if (!mumps->sym) { /* MUMPS always return a full matrix */ 18166444a565SStefano Zampini if (mumps->id.ICNTL(19) == 1) { /* stored by rows */ 18176444a565SStefano Zampini PetscInt i,j,N=mumps->id.size_schur; 18186444a565SStefano Zampini for (i=0;i<N;i++) { 18196444a565SStefano Zampini for (j=0;j<N;j++) { 18206444a565SStefano Zampini #if !defined(PETSC_USE_COMPLEX) 18216444a565SStefano Zampini PetscScalar val = mumps->id.schur[i*N+j]; 18226444a565SStefano Zampini #else 18236444a565SStefano Zampini PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i; 18246444a565SStefano Zampini #endif 18256444a565SStefano Zampini array[j*N+i] = val; 18266444a565SStefano Zampini } 18276444a565SStefano Zampini } 18286444a565SStefano Zampini } else { /* stored by columns */ 18296444a565SStefano Zampini ierr = PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));CHKERRQ(ierr); 18306444a565SStefano Zampini } 18316444a565SStefano Zampini } else { /* either full or lower-triangular (not packed) */ 18326444a565SStefano Zampini if (mumps->id.ICNTL(19) == 2) { /* lower triangular stored by columns */ 18336444a565SStefano Zampini PetscInt i,j,N=mumps->id.size_schur; 18346444a565SStefano Zampini for (i=0;i<N;i++) { 18356444a565SStefano Zampini for (j=i;j<N;j++) { 18366444a565SStefano Zampini #if !defined(PETSC_USE_COMPLEX) 18376444a565SStefano Zampini PetscScalar val = mumps->id.schur[i*N+j]; 18386444a565SStefano Zampini #else 18396444a565SStefano Zampini PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i; 18406444a565SStefano Zampini #endif 18416444a565SStefano Zampini array[i*N+j] = val; 18426444a565SStefano Zampini array[j*N+i] = val; 18436444a565SStefano Zampini } 18446444a565SStefano Zampini } 18456444a565SStefano Zampini } else if (mumps->id.ICNTL(19) == 3) { /* full matrix */ 18466444a565SStefano Zampini ierr = PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));CHKERRQ(ierr); 18476444a565SStefano Zampini } else { /* ICNTL(19) == 1 lower triangular stored by rows */ 18486444a565SStefano Zampini PetscInt i,j,N=mumps->id.size_schur; 18496444a565SStefano Zampini for (i=0;i<N;i++) { 18506444a565SStefano Zampini for (j=0;j<i+1;j++) { 18516444a565SStefano Zampini #if !defined(PETSC_USE_COMPLEX) 18526444a565SStefano Zampini PetscScalar val = mumps->id.schur[i*N+j]; 18536444a565SStefano Zampini #else 18546444a565SStefano Zampini PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i; 18556444a565SStefano Zampini #endif 18566444a565SStefano Zampini array[i*N+j] = val; 18576444a565SStefano Zampini array[j*N+i] = val; 18586444a565SStefano Zampini } 18596444a565SStefano Zampini } 18606444a565SStefano Zampini } 18616444a565SStefano Zampini } 18626444a565SStefano Zampini ierr = MatDenseRestoreArray(St,&array);CHKERRQ(ierr); 18636444a565SStefano Zampini *S = St; 18646444a565SStefano Zampini PetscFunctionReturn(0); 18656444a565SStefano Zampini } 18666444a565SStefano Zampini 186759ac8732SStefano Zampini /* -------------------------------------------------------------------------------------------*/ 186859ac8732SStefano Zampini #undef __FUNCT__ 18695a05ddb0SStefano Zampini #define __FUNCT__ "MatFactorGetSchurComplement_MUMPS" 18705a05ddb0SStefano Zampini PetscErrorCode MatFactorGetSchurComplement_MUMPS(Mat F,Mat* S) 187159ac8732SStefano Zampini { 187259ac8732SStefano Zampini Mat St; 187359ac8732SStefano Zampini Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 187459ac8732SStefano Zampini PetscErrorCode ierr; 187559ac8732SStefano Zampini 187659ac8732SStefano Zampini PetscFunctionBegin; 18775a05ddb0SStefano 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"); 18785a05ddb0SStefano Zampini else if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before"); 18796f3cc6f9SBarry Smith 18807bc32933SStefano 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 */ 188159ac8732SStefano Zampini ierr = MatCreateSeqDense(PetscObjectComm((PetscObject)F),mumps->id.size_schur,mumps->id.size_schur,(PetscScalar*)mumps->id.schur,&St);CHKERRQ(ierr); 188259ac8732SStefano Zampini *S = St; 188359ac8732SStefano Zampini PetscFunctionReturn(0); 188459ac8732SStefano Zampini } 188559ac8732SStefano Zampini 188659ac8732SStefano Zampini /* -------------------------------------------------------------------------------------------*/ 188759ac8732SStefano Zampini #undef __FUNCT__ 18885a05ddb0SStefano Zampini #define __FUNCT__ "MatFactorInvertSchurComplement_MUMPS" 18895a05ddb0SStefano Zampini PetscErrorCode MatFactorInvertSchurComplement_MUMPS(Mat F) 189059ac8732SStefano Zampini { 189159ac8732SStefano Zampini Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 189259ac8732SStefano Zampini PetscErrorCode ierr; 189359ac8732SStefano Zampini 189459ac8732SStefano Zampini PetscFunctionBegin; 189559ac8732SStefano Zampini if (!mumps->id.ICNTL(19)) { /* do nothing */ 189659ac8732SStefano Zampini PetscFunctionReturn(0); 189759ac8732SStefano Zampini } 18985a05ddb0SStefano Zampini if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before"); 189959ac8732SStefano Zampini ierr = MatMumpsInvertSchur_Private(mumps);CHKERRQ(ierr); 190059ac8732SStefano Zampini PetscFunctionReturn(0); 190159ac8732SStefano Zampini } 190259ac8732SStefano Zampini 19036444a565SStefano Zampini /* -------------------------------------------------------------------------------------------*/ 19046444a565SStefano Zampini #undef __FUNCT__ 19055a05ddb0SStefano Zampini #define __FUNCT__ "MatFactorSolveSchurComplement_MUMPS" 19065a05ddb0SStefano Zampini PetscErrorCode MatFactorSolveSchurComplement_MUMPS(Mat F, Vec rhs, Vec sol) 1907e807eca7SStefano Zampini { 1908e807eca7SStefano Zampini Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 1909e807eca7SStefano Zampini MumpsScalar *orhs; 1910e807eca7SStefano Zampini PetscScalar *osol,*nrhs,*nsol; 19119a3a5937SStefano Zampini PetscInt orhs_size,osol_size,olrhs_size; 1912e807eca7SStefano Zampini PetscErrorCode ierr; 1913e807eca7SStefano Zampini 1914e807eca7SStefano Zampini PetscFunctionBegin; 19155a05ddb0SStefano 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"); 19165a05ddb0SStefano Zampini if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before"); 19176f3cc6f9SBarry Smith 1918e807eca7SStefano Zampini /* swap pointers */ 1919e807eca7SStefano Zampini orhs = mumps->id.redrhs; 19209a3a5937SStefano Zampini olrhs_size = mumps->id.lredrhs; 1921a12f35bfSStefano Zampini orhs_size = mumps->sizeredrhs; 1922e807eca7SStefano Zampini osol = mumps->schur_sol; 1923a12f35bfSStefano Zampini osol_size = mumps->schur_sizesol; 1924e807eca7SStefano Zampini ierr = VecGetArray(rhs,&nrhs);CHKERRQ(ierr); 1925e807eca7SStefano Zampini ierr = VecGetArray(sol,&nsol);CHKERRQ(ierr); 1926e807eca7SStefano Zampini mumps->id.redrhs = (MumpsScalar*)nrhs; 1927a12f35bfSStefano Zampini ierr = VecGetLocalSize(rhs,&mumps->sizeredrhs);CHKERRQ(ierr); 19289a3a5937SStefano Zampini mumps->id.lredrhs = mumps->sizeredrhs; 1929e807eca7SStefano Zampini mumps->schur_sol = nsol; 1930a12f35bfSStefano Zampini ierr = VecGetLocalSize(sol,&mumps->schur_sizesol);CHKERRQ(ierr); 1931a12f35bfSStefano Zampini 1932e807eca7SStefano Zampini /* solve Schur complement */ 1933e807eca7SStefano Zampini mumps->id.nrhs = 1; 1934e807eca7SStefano Zampini ierr = MatMumpsSolveSchur_Private(mumps,PETSC_FALSE);CHKERRQ(ierr); 1935e807eca7SStefano Zampini /* restore pointers */ 1936e807eca7SStefano Zampini ierr = VecRestoreArray(rhs,&nrhs);CHKERRQ(ierr); 1937e807eca7SStefano Zampini ierr = VecRestoreArray(sol,&nsol);CHKERRQ(ierr); 1938e807eca7SStefano Zampini mumps->id.redrhs = orhs; 19399a3a5937SStefano Zampini mumps->id.lredrhs = olrhs_size; 1940a12f35bfSStefano Zampini mumps->sizeredrhs = orhs_size; 1941e807eca7SStefano Zampini mumps->schur_sol = osol; 1942a12f35bfSStefano Zampini mumps->schur_sizesol = osol_size; 1943e807eca7SStefano Zampini PetscFunctionReturn(0); 1944e807eca7SStefano Zampini } 1945e807eca7SStefano Zampini 19467404bcfbSStefano Zampini /* -------------------------------------------------------------------------------------------*/ 19477404bcfbSStefano Zampini #undef __FUNCT__ 19485a05ddb0SStefano Zampini #define __FUNCT__ "MatFactorSolveSchurComplementTranspose_MUMPS" 19495a05ddb0SStefano Zampini PetscErrorCode MatFactorSolveSchurComplementTranspose_MUMPS(Mat F, Vec rhs, Vec sol) 19507404bcfbSStefano Zampini { 19517404bcfbSStefano Zampini Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 19527404bcfbSStefano Zampini MumpsScalar *orhs; 19537404bcfbSStefano Zampini PetscScalar *osol,*nrhs,*nsol; 1954a12f35bfSStefano Zampini PetscInt orhs_size,osol_size; 19557404bcfbSStefano Zampini PetscErrorCode ierr; 19567404bcfbSStefano Zampini 19577404bcfbSStefano Zampini PetscFunctionBegin; 19585a05ddb0SStefano 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"); 19595a05ddb0SStefano Zampini else if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before"); 19606f3cc6f9SBarry Smith 19617404bcfbSStefano Zampini /* swap pointers */ 19627404bcfbSStefano Zampini orhs = mumps->id.redrhs; 1963a12f35bfSStefano Zampini orhs_size = mumps->sizeredrhs; 19647404bcfbSStefano Zampini osol = mumps->schur_sol; 1965a12f35bfSStefano Zampini osol_size = mumps->schur_sizesol; 19667404bcfbSStefano Zampini ierr = VecGetArray(rhs,&nrhs);CHKERRQ(ierr); 19677404bcfbSStefano Zampini ierr = VecGetArray(sol,&nsol);CHKERRQ(ierr); 19687404bcfbSStefano Zampini mumps->id.redrhs = (MumpsScalar*)nrhs; 1969a12f35bfSStefano Zampini ierr = VecGetLocalSize(rhs,&mumps->sizeredrhs);CHKERRQ(ierr); 19707404bcfbSStefano Zampini mumps->schur_sol = nsol; 1971a12f35bfSStefano Zampini ierr = VecGetLocalSize(sol,&mumps->schur_sizesol);CHKERRQ(ierr); 1972a12f35bfSStefano Zampini 19737404bcfbSStefano Zampini /* solve Schur complement */ 19747404bcfbSStefano Zampini mumps->id.nrhs = 1; 19757404bcfbSStefano Zampini mumps->id.ICNTL(9) = 0; 19767404bcfbSStefano Zampini ierr = MatMumpsSolveSchur_Private(mumps,PETSC_FALSE);CHKERRQ(ierr); 19777404bcfbSStefano Zampini mumps->id.ICNTL(9) = 1; 19787404bcfbSStefano Zampini /* restore pointers */ 19797404bcfbSStefano Zampini ierr = VecRestoreArray(rhs,&nrhs);CHKERRQ(ierr); 19807404bcfbSStefano Zampini ierr = VecRestoreArray(sol,&nsol);CHKERRQ(ierr); 19817404bcfbSStefano Zampini mumps->id.redrhs = orhs; 1982a12f35bfSStefano Zampini mumps->sizeredrhs = orhs_size; 19837404bcfbSStefano Zampini mumps->schur_sol = osol; 1984a12f35bfSStefano Zampini mumps->schur_sizesol = osol_size; 19857404bcfbSStefano Zampini PetscFunctionReturn(0); 19867404bcfbSStefano Zampini } 19877404bcfbSStefano Zampini 1988e807eca7SStefano Zampini /* -------------------------------------------------------------------------------------------*/ 1989e807eca7SStefano Zampini #undef __FUNCT__ 19905ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl_MUMPS" 19915ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival) 19925ccb76cbSHong Zhang { 1993a5e57a09SHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 19945ccb76cbSHong Zhang 19955ccb76cbSHong Zhang PetscFunctionBegin; 1996a5e57a09SHong Zhang mumps->id.ICNTL(icntl) = ival; 19975ccb76cbSHong Zhang PetscFunctionReturn(0); 19985ccb76cbSHong Zhang } 19995ccb76cbSHong Zhang 20005ccb76cbSHong Zhang #undef __FUNCT__ 2001bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetIcntl_MUMPS" 2002bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt *ival) 2003bc6112feSHong Zhang { 2004bc6112feSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 2005bc6112feSHong Zhang 2006bc6112feSHong Zhang PetscFunctionBegin; 2007bc6112feSHong Zhang *ival = mumps->id.ICNTL(icntl); 2008bc6112feSHong Zhang PetscFunctionReturn(0); 2009bc6112feSHong Zhang } 2010bc6112feSHong Zhang 2011bc6112feSHong Zhang #undef __FUNCT__ 20125ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl" 20135ccb76cbSHong Zhang /*@ 20145ccb76cbSHong Zhang MatMumpsSetIcntl - Set MUMPS parameter ICNTL() 20155ccb76cbSHong Zhang 20165ccb76cbSHong Zhang Logically Collective on Mat 20175ccb76cbSHong Zhang 20185ccb76cbSHong Zhang Input Parameters: 20195ccb76cbSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 20205ccb76cbSHong Zhang . icntl - index of MUMPS parameter array ICNTL() 20215ccb76cbSHong Zhang - ival - value of MUMPS ICNTL(icntl) 20225ccb76cbSHong Zhang 20235ccb76cbSHong Zhang Options Database: 20245ccb76cbSHong Zhang . -mat_mumps_icntl_<icntl> <ival> 20255ccb76cbSHong Zhang 20265ccb76cbSHong Zhang Level: beginner 20275ccb76cbSHong Zhang 202896a0c994SBarry Smith References: 202996a0c994SBarry Smith . MUMPS Users' Guide 20305ccb76cbSHong Zhang 20315ccb76cbSHong Zhang .seealso: MatGetFactor() 20325ccb76cbSHong Zhang @*/ 20335ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival) 20345ccb76cbSHong Zhang { 20355ccb76cbSHong Zhang PetscErrorCode ierr; 20365ccb76cbSHong Zhang 20375ccb76cbSHong Zhang PetscFunctionBegin; 2038*2989dfd4SHong Zhang PetscValidType(F,1); 2039*2989dfd4SHong Zhang if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix"); 20405ccb76cbSHong Zhang PetscValidLogicalCollectiveInt(F,icntl,2); 20415ccb76cbSHong Zhang PetscValidLogicalCollectiveInt(F,ival,3); 20425ccb76cbSHong Zhang ierr = PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr); 20435ccb76cbSHong Zhang PetscFunctionReturn(0); 20445ccb76cbSHong Zhang } 20455ccb76cbSHong Zhang 2046bc6112feSHong Zhang #undef __FUNCT__ 2047bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetIcntl" 2048a21f80fcSHong Zhang /*@ 2049a21f80fcSHong Zhang MatMumpsGetIcntl - Get MUMPS parameter ICNTL() 2050a21f80fcSHong Zhang 2051a21f80fcSHong Zhang Logically Collective on Mat 2052a21f80fcSHong Zhang 2053a21f80fcSHong Zhang Input Parameters: 2054a21f80fcSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 2055a21f80fcSHong Zhang - icntl - index of MUMPS parameter array ICNTL() 2056a21f80fcSHong Zhang 2057a21f80fcSHong Zhang Output Parameter: 2058a21f80fcSHong Zhang . ival - value of MUMPS ICNTL(icntl) 2059a21f80fcSHong Zhang 2060a21f80fcSHong Zhang Level: beginner 2061a21f80fcSHong Zhang 206296a0c994SBarry Smith References: 206396a0c994SBarry Smith . MUMPS Users' Guide 2064a21f80fcSHong Zhang 2065a21f80fcSHong Zhang .seealso: MatGetFactor() 2066a21f80fcSHong Zhang @*/ 2067bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl(Mat F,PetscInt icntl,PetscInt *ival) 2068bc6112feSHong Zhang { 2069bc6112feSHong Zhang PetscErrorCode ierr; 2070bc6112feSHong Zhang 2071bc6112feSHong Zhang PetscFunctionBegin; 2072*2989dfd4SHong Zhang PetscValidType(F,1); 2073*2989dfd4SHong Zhang if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix"); 2074bc6112feSHong Zhang PetscValidLogicalCollectiveInt(F,icntl,2); 2075bc6112feSHong Zhang PetscValidIntPointer(ival,3); 2076*2989dfd4SHong Zhang ierr = PetscUseMethod(F,"MatMumpsGetIcntl_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr); 2077bc6112feSHong Zhang PetscFunctionReturn(0); 2078bc6112feSHong Zhang } 2079bc6112feSHong Zhang 20808928b65cSHong Zhang /* -------------------------------------------------------------------------------------------*/ 20818928b65cSHong Zhang #undef __FUNCT__ 20828928b65cSHong Zhang #define __FUNCT__ "MatMumpsSetCntl_MUMPS" 20838928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal val) 20848928b65cSHong Zhang { 20858928b65cSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 20868928b65cSHong Zhang 20878928b65cSHong Zhang PetscFunctionBegin; 20888928b65cSHong Zhang mumps->id.CNTL(icntl) = val; 20898928b65cSHong Zhang PetscFunctionReturn(0); 20908928b65cSHong Zhang } 20918928b65cSHong Zhang 20928928b65cSHong Zhang #undef __FUNCT__ 2093bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetCntl_MUMPS" 2094bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal *val) 2095bc6112feSHong Zhang { 2096bc6112feSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 2097bc6112feSHong Zhang 2098bc6112feSHong Zhang PetscFunctionBegin; 2099bc6112feSHong Zhang *val = mumps->id.CNTL(icntl); 2100bc6112feSHong Zhang PetscFunctionReturn(0); 2101bc6112feSHong Zhang } 2102bc6112feSHong Zhang 2103bc6112feSHong Zhang #undef __FUNCT__ 21048928b65cSHong Zhang #define __FUNCT__ "MatMumpsSetCntl" 21058928b65cSHong Zhang /*@ 21068928b65cSHong Zhang MatMumpsSetCntl - Set MUMPS parameter CNTL() 21078928b65cSHong Zhang 21088928b65cSHong Zhang Logically Collective on Mat 21098928b65cSHong Zhang 21108928b65cSHong Zhang Input Parameters: 21118928b65cSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 21128928b65cSHong Zhang . icntl - index of MUMPS parameter array CNTL() 21138928b65cSHong Zhang - val - value of MUMPS CNTL(icntl) 21148928b65cSHong Zhang 21158928b65cSHong Zhang Options Database: 21168928b65cSHong Zhang . -mat_mumps_cntl_<icntl> <val> 21178928b65cSHong Zhang 21188928b65cSHong Zhang Level: beginner 21198928b65cSHong Zhang 212096a0c994SBarry Smith References: 212196a0c994SBarry Smith . MUMPS Users' Guide 21228928b65cSHong Zhang 21238928b65cSHong Zhang .seealso: MatGetFactor() 21248928b65cSHong Zhang @*/ 21258928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl(Mat F,PetscInt icntl,PetscReal val) 21268928b65cSHong Zhang { 21278928b65cSHong Zhang PetscErrorCode ierr; 21288928b65cSHong Zhang 21298928b65cSHong Zhang PetscFunctionBegin; 2130*2989dfd4SHong Zhang PetscValidType(F,1); 2131*2989dfd4SHong Zhang if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix"); 21328928b65cSHong Zhang PetscValidLogicalCollectiveInt(F,icntl,2); 2133bc6112feSHong Zhang PetscValidLogicalCollectiveReal(F,val,3); 21348928b65cSHong Zhang ierr = PetscTryMethod(F,"MatMumpsSetCntl_C",(Mat,PetscInt,PetscReal),(F,icntl,val));CHKERRQ(ierr); 21358928b65cSHong Zhang PetscFunctionReturn(0); 21368928b65cSHong Zhang } 21378928b65cSHong Zhang 2138bc6112feSHong Zhang #undef __FUNCT__ 2139bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetCntl" 2140a21f80fcSHong Zhang /*@ 2141a21f80fcSHong Zhang MatMumpsGetCntl - Get MUMPS parameter CNTL() 2142a21f80fcSHong Zhang 2143a21f80fcSHong Zhang Logically Collective on Mat 2144a21f80fcSHong Zhang 2145a21f80fcSHong Zhang Input Parameters: 2146a21f80fcSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 2147a21f80fcSHong Zhang - icntl - index of MUMPS parameter array CNTL() 2148a21f80fcSHong Zhang 2149a21f80fcSHong Zhang Output Parameter: 2150a21f80fcSHong Zhang . val - value of MUMPS CNTL(icntl) 2151a21f80fcSHong Zhang 2152a21f80fcSHong Zhang Level: beginner 2153a21f80fcSHong Zhang 215496a0c994SBarry Smith References: 215596a0c994SBarry Smith . MUMPS Users' Guide 2156a21f80fcSHong Zhang 2157a21f80fcSHong Zhang .seealso: MatGetFactor() 2158a21f80fcSHong Zhang @*/ 2159bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl(Mat F,PetscInt icntl,PetscReal *val) 2160bc6112feSHong Zhang { 2161bc6112feSHong Zhang PetscErrorCode ierr; 2162bc6112feSHong Zhang 2163bc6112feSHong Zhang PetscFunctionBegin; 2164*2989dfd4SHong Zhang PetscValidType(F,1); 2165*2989dfd4SHong Zhang if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix"); 2166bc6112feSHong Zhang PetscValidLogicalCollectiveInt(F,icntl,2); 2167bc6112feSHong Zhang PetscValidRealPointer(val,3); 2168*2989dfd4SHong Zhang ierr = PetscUseMethod(F,"MatMumpsGetCntl_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr); 2169bc6112feSHong Zhang PetscFunctionReturn(0); 2170bc6112feSHong Zhang } 2171bc6112feSHong Zhang 2172bc6112feSHong Zhang #undef __FUNCT__ 2173ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfo_MUMPS" 2174ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F,PetscInt icntl,PetscInt *info) 2175bc6112feSHong Zhang { 2176bc6112feSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 2177bc6112feSHong Zhang 2178bc6112feSHong Zhang PetscFunctionBegin; 2179bc6112feSHong Zhang *info = mumps->id.INFO(icntl); 2180bc6112feSHong Zhang PetscFunctionReturn(0); 2181bc6112feSHong Zhang } 2182bc6112feSHong Zhang 2183bc6112feSHong Zhang #undef __FUNCT__ 2184ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfog_MUMPS" 2185ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F,PetscInt icntl,PetscInt *infog) 2186bc6112feSHong Zhang { 2187bc6112feSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 2188bc6112feSHong Zhang 2189bc6112feSHong Zhang PetscFunctionBegin; 2190bc6112feSHong Zhang *infog = mumps->id.INFOG(icntl); 2191bc6112feSHong Zhang PetscFunctionReturn(0); 2192bc6112feSHong Zhang } 2193bc6112feSHong Zhang 2194bc6112feSHong Zhang #undef __FUNCT__ 2195ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfo_MUMPS" 2196ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfo) 2197bc6112feSHong Zhang { 2198bc6112feSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 2199bc6112feSHong Zhang 2200bc6112feSHong Zhang PetscFunctionBegin; 2201bc6112feSHong Zhang *rinfo = mumps->id.RINFO(icntl); 2202bc6112feSHong Zhang PetscFunctionReturn(0); 2203bc6112feSHong Zhang } 2204bc6112feSHong Zhang 2205bc6112feSHong Zhang #undef __FUNCT__ 2206ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfog_MUMPS" 2207ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfog) 2208bc6112feSHong Zhang { 2209bc6112feSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 2210bc6112feSHong Zhang 2211bc6112feSHong Zhang PetscFunctionBegin; 2212bc6112feSHong Zhang *rinfog = mumps->id.RINFOG(icntl); 2213bc6112feSHong Zhang PetscFunctionReturn(0); 2214bc6112feSHong Zhang } 2215bc6112feSHong Zhang 2216bc6112feSHong Zhang #undef __FUNCT__ 2217ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfo" 2218a21f80fcSHong Zhang /*@ 2219a21f80fcSHong Zhang MatMumpsGetInfo - Get MUMPS parameter INFO() 2220a21f80fcSHong Zhang 2221a21f80fcSHong Zhang Logically Collective on Mat 2222a21f80fcSHong Zhang 2223a21f80fcSHong Zhang Input Parameters: 2224a21f80fcSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 2225a21f80fcSHong Zhang - icntl - index of MUMPS parameter array INFO() 2226a21f80fcSHong Zhang 2227a21f80fcSHong Zhang Output Parameter: 2228a21f80fcSHong Zhang . ival - value of MUMPS INFO(icntl) 2229a21f80fcSHong Zhang 2230a21f80fcSHong Zhang Level: beginner 2231a21f80fcSHong Zhang 223296a0c994SBarry Smith References: 223396a0c994SBarry Smith . MUMPS Users' Guide 2234a21f80fcSHong Zhang 2235a21f80fcSHong Zhang .seealso: MatGetFactor() 2236a21f80fcSHong Zhang @*/ 2237ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo(Mat F,PetscInt icntl,PetscInt *ival) 2238bc6112feSHong Zhang { 2239bc6112feSHong Zhang PetscErrorCode ierr; 2240bc6112feSHong Zhang 2241bc6112feSHong Zhang PetscFunctionBegin; 2242*2989dfd4SHong Zhang PetscValidType(F,1); 2243*2989dfd4SHong Zhang if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix"); 2244ca810319SHong Zhang PetscValidIntPointer(ival,3); 2245*2989dfd4SHong Zhang ierr = PetscUseMethod(F,"MatMumpsGetInfo_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr); 2246bc6112feSHong Zhang PetscFunctionReturn(0); 2247bc6112feSHong Zhang } 2248bc6112feSHong Zhang 2249bc6112feSHong Zhang #undef __FUNCT__ 2250ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfog" 2251a21f80fcSHong Zhang /*@ 2252a21f80fcSHong Zhang MatMumpsGetInfog - Get MUMPS parameter INFOG() 2253a21f80fcSHong Zhang 2254a21f80fcSHong Zhang Logically Collective on Mat 2255a21f80fcSHong Zhang 2256a21f80fcSHong Zhang Input Parameters: 2257a21f80fcSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 2258a21f80fcSHong Zhang - icntl - index of MUMPS parameter array INFOG() 2259a21f80fcSHong Zhang 2260a21f80fcSHong Zhang Output Parameter: 2261a21f80fcSHong Zhang . ival - value of MUMPS INFOG(icntl) 2262a21f80fcSHong Zhang 2263a21f80fcSHong Zhang Level: beginner 2264a21f80fcSHong Zhang 226596a0c994SBarry Smith References: 226696a0c994SBarry Smith . MUMPS Users' Guide 2267a21f80fcSHong Zhang 2268a21f80fcSHong Zhang .seealso: MatGetFactor() 2269a21f80fcSHong Zhang @*/ 2270ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog(Mat F,PetscInt icntl,PetscInt *ival) 2271bc6112feSHong Zhang { 2272bc6112feSHong Zhang PetscErrorCode ierr; 2273bc6112feSHong Zhang 2274bc6112feSHong Zhang PetscFunctionBegin; 2275*2989dfd4SHong Zhang PetscValidType(F,1); 2276*2989dfd4SHong Zhang if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix"); 2277ca810319SHong Zhang PetscValidIntPointer(ival,3); 2278*2989dfd4SHong Zhang ierr = PetscUseMethod(F,"MatMumpsGetInfog_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr); 2279bc6112feSHong Zhang PetscFunctionReturn(0); 2280bc6112feSHong Zhang } 2281bc6112feSHong Zhang 2282bc6112feSHong Zhang #undef __FUNCT__ 2283ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfo" 2284a21f80fcSHong Zhang /*@ 2285a21f80fcSHong Zhang MatMumpsGetRinfo - Get MUMPS parameter RINFO() 2286a21f80fcSHong Zhang 2287a21f80fcSHong Zhang Logically Collective on Mat 2288a21f80fcSHong Zhang 2289a21f80fcSHong Zhang Input Parameters: 2290a21f80fcSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 2291a21f80fcSHong Zhang - icntl - index of MUMPS parameter array RINFO() 2292a21f80fcSHong Zhang 2293a21f80fcSHong Zhang Output Parameter: 2294a21f80fcSHong Zhang . val - value of MUMPS RINFO(icntl) 2295a21f80fcSHong Zhang 2296a21f80fcSHong Zhang Level: beginner 2297a21f80fcSHong Zhang 229896a0c994SBarry Smith References: 229996a0c994SBarry Smith . MUMPS Users' Guide 2300a21f80fcSHong Zhang 2301a21f80fcSHong Zhang .seealso: MatGetFactor() 2302a21f80fcSHong Zhang @*/ 2303ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo(Mat F,PetscInt icntl,PetscReal *val) 2304bc6112feSHong Zhang { 2305bc6112feSHong Zhang PetscErrorCode ierr; 2306bc6112feSHong Zhang 2307bc6112feSHong Zhang PetscFunctionBegin; 2308*2989dfd4SHong Zhang PetscValidType(F,1); 2309*2989dfd4SHong Zhang if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix"); 2310bc6112feSHong Zhang PetscValidRealPointer(val,3); 2311*2989dfd4SHong Zhang ierr = PetscUseMethod(F,"MatMumpsGetRinfo_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr); 2312bc6112feSHong Zhang PetscFunctionReturn(0); 2313bc6112feSHong Zhang } 2314bc6112feSHong Zhang 2315bc6112feSHong Zhang #undef __FUNCT__ 2316ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfog" 2317a21f80fcSHong Zhang /*@ 2318a21f80fcSHong Zhang MatMumpsGetRinfog - Get MUMPS parameter RINFOG() 2319a21f80fcSHong Zhang 2320a21f80fcSHong Zhang Logically Collective on Mat 2321a21f80fcSHong Zhang 2322a21f80fcSHong Zhang Input Parameters: 2323a21f80fcSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 2324a21f80fcSHong Zhang - icntl - index of MUMPS parameter array RINFOG() 2325a21f80fcSHong Zhang 2326a21f80fcSHong Zhang Output Parameter: 2327a21f80fcSHong Zhang . val - value of MUMPS RINFOG(icntl) 2328a21f80fcSHong Zhang 2329a21f80fcSHong Zhang Level: beginner 2330a21f80fcSHong Zhang 233196a0c994SBarry Smith References: 233296a0c994SBarry Smith . MUMPS Users' Guide 2333a21f80fcSHong Zhang 2334a21f80fcSHong Zhang .seealso: MatGetFactor() 2335a21f80fcSHong Zhang @*/ 2336ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog(Mat F,PetscInt icntl,PetscReal *val) 2337bc6112feSHong Zhang { 2338bc6112feSHong Zhang PetscErrorCode ierr; 2339bc6112feSHong Zhang 2340bc6112feSHong Zhang PetscFunctionBegin; 2341*2989dfd4SHong Zhang PetscValidType(F,1); 2342*2989dfd4SHong Zhang if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix"); 2343bc6112feSHong Zhang PetscValidRealPointer(val,3); 2344*2989dfd4SHong Zhang ierr = PetscUseMethod(F,"MatMumpsGetRinfog_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr); 2345bc6112feSHong Zhang PetscFunctionReturn(0); 2346bc6112feSHong Zhang } 2347bc6112feSHong Zhang 234824b6179bSKris Buschelman /*MC 23492692d6eeSBarry Smith MATSOLVERMUMPS - A matrix type providing direct solvers (LU and Cholesky) for 235024b6179bSKris Buschelman distributed and sequential matrices via the external package MUMPS. 235124b6179bSKris Buschelman 235241c8de11SBarry Smith Works with MATAIJ and MATSBAIJ matrices 235324b6179bSKris Buschelman 2354c2b89b5dSBarry Smith Use ./configure --download-mumps --download-scalapack --download-parmetis --download-metis --download-ptscotch to have PETSc installed with MUMPS 2355c2b89b5dSBarry Smith 2356c2b89b5dSBarry Smith Use -pc_type cholesky or lu -pc_factor_mat_solver_package mumps to us this direct solver 2357c2b89b5dSBarry Smith 235824b6179bSKris Buschelman Options Database Keys: 23594e34a73bSHong Zhang + -mat_mumps_icntl_1 <6>: ICNTL(1): output stream for error messages (None) 23604e34a73bSHong Zhang . -mat_mumps_icntl_2 <0>: ICNTL(2): output stream for diagnostic printing, statistics, and warning (None) 23614e34a73bSHong Zhang . -mat_mumps_icntl_3 <0>: ICNTL(3): output stream for global information, collected on the host (None) 23624e34a73bSHong Zhang . -mat_mumps_icntl_4 <0>: ICNTL(4): level of printing (0 to 4) (None) 23634e34a73bSHong Zhang . -mat_mumps_icntl_6 <7>: ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7) (None) 23644e34a73bSHong 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) 23654e34a73bSHong Zhang . -mat_mumps_icntl_8 <77>: ICNTL(8): scaling strategy (-2 to 8 or 77) (None) 23664e34a73bSHong Zhang . -mat_mumps_icntl_10 <0>: ICNTL(10): max num of refinements (None) 23674e34a73bSHong Zhang . -mat_mumps_icntl_11 <0>: ICNTL(11): statistics related to an error analysis (via -ksp_view) (None) 23684e34a73bSHong Zhang . -mat_mumps_icntl_12 <1>: ICNTL(12): an ordering strategy for symmetric matrices (0 to 3) (None) 23694e34a73bSHong Zhang . -mat_mumps_icntl_13 <0>: ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting (None) 23704e34a73bSHong Zhang . -mat_mumps_icntl_14 <20>: ICNTL(14): percentage increase in the estimated working space (None) 23714e34a73bSHong Zhang . -mat_mumps_icntl_19 <0>: ICNTL(19): computes the Schur complement (None) 23724e34a73bSHong Zhang . -mat_mumps_icntl_22 <0>: ICNTL(22): in-core/out-of-core factorization and solve (0 or 1) (None) 23734e34a73bSHong Zhang . -mat_mumps_icntl_23 <0>: ICNTL(23): max size of the working memory (MB) that can allocate per processor (None) 23744e34a73bSHong Zhang . -mat_mumps_icntl_24 <0>: ICNTL(24): detection of null pivot rows (0 or 1) (None) 23754e34a73bSHong Zhang . -mat_mumps_icntl_25 <0>: ICNTL(25): compute a solution of a deficient matrix and a null space basis (None) 23764e34a73bSHong Zhang . -mat_mumps_icntl_26 <0>: ICNTL(26): drives the solution phase if a Schur complement matrix (None) 23774e34a73bSHong 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) 23784e34a73bSHong Zhang . -mat_mumps_icntl_29 <0>: ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis (None) 23794e34a73bSHong Zhang . -mat_mumps_icntl_30 <0>: ICNTL(30): compute user-specified set of entries in inv(A) (None) 23804e34a73bSHong Zhang . -mat_mumps_icntl_31 <0>: ICNTL(31): indicates which factors may be discarded during factorization (None) 23814e34a73bSHong Zhang . -mat_mumps_icntl_33 <0>: ICNTL(33): compute determinant (None) 23824e34a73bSHong Zhang . -mat_mumps_cntl_1 <0.01>: CNTL(1): relative pivoting threshold (None) 23834e34a73bSHong Zhang . -mat_mumps_cntl_2 <1.49012e-08>: CNTL(2): stopping criterion of refinement (None) 23844e34a73bSHong Zhang . -mat_mumps_cntl_3 <0>: CNTL(3): absolute pivoting threshold (None) 23854e34a73bSHong Zhang . -mat_mumps_cntl_4 <-1>: CNTL(4): value for static pivoting (None) 23864e34a73bSHong Zhang - -mat_mumps_cntl_5 <0>: CNTL(5): fixation for null pivots (None) 238724b6179bSKris Buschelman 238824b6179bSKris Buschelman Level: beginner 238924b6179bSKris Buschelman 239041c8de11SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage 239141c8de11SBarry Smith 239224b6179bSKris Buschelman M*/ 239324b6179bSKris Buschelman 239435bd34faSBarry Smith #undef __FUNCT__ 239535bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_mumps" 2396f7a08781SBarry Smith static PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type) 239735bd34faSBarry Smith { 239835bd34faSBarry Smith PetscFunctionBegin; 23992692d6eeSBarry Smith *type = MATSOLVERMUMPS; 240035bd34faSBarry Smith PetscFunctionReturn(0); 240135bd34faSBarry Smith } 240235bd34faSBarry Smith 2403bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI AIJ matrices */ 24042877fffaSHong Zhang #undef __FUNCT__ 2405bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_aij_mumps" 24068cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F) 24072877fffaSHong Zhang { 24082877fffaSHong Zhang Mat B; 24092877fffaSHong Zhang PetscErrorCode ierr; 24102877fffaSHong Zhang Mat_MUMPS *mumps; 2411ace3abfcSBarry Smith PetscBool isSeqAIJ; 24122877fffaSHong Zhang 24132877fffaSHong Zhang PetscFunctionBegin; 24142877fffaSHong Zhang /* Create the factorization matrix */ 2415251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr); 2416ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 24172877fffaSHong Zhang ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 24182877fffaSHong Zhang ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 2419bccb9932SShri Abhyankar if (isSeqAIJ) { 24200298fd71SBarry Smith ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr); 2421bccb9932SShri Abhyankar } else { 24220298fd71SBarry Smith ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr); 2423bccb9932SShri Abhyankar } 24242877fffaSHong Zhang 2425b00a9115SJed Brown ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr); 24262205254eSKarl Rupp 24272877fffaSHong Zhang B->ops->view = MatView_MUMPS; 242835bd34faSBarry Smith B->ops->getinfo = MatGetInfo_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); 2438bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr); 2439bc6112feSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr); 2440bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr); 2441bc6112feSHong 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); 24466444a565SStefano Zampini 2447450b117fSShri Abhyankar if (ftype == MAT_FACTOR_LU) { 2448450b117fSShri Abhyankar B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 2449d5f3da31SBarry Smith B->factortype = MAT_FACTOR_LU; 2450bccb9932SShri Abhyankar if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij; 2451bccb9932SShri Abhyankar else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij; 2452746480a1SHong Zhang mumps->sym = 0; 2453dcd589f8SShri Abhyankar } else { 245467877ebaSShri Abhyankar B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS; 2455450b117fSShri Abhyankar B->factortype = MAT_FACTOR_CHOLESKY; 2456bccb9932SShri Abhyankar if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij; 2457bccb9932SShri Abhyankar else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij; 245859ac8732SStefano Zampini #if defined(PETSC_USE_COMPLEX) 245959ac8732SStefano Zampini mumps->sym = 2; 246059ac8732SStefano Zampini #else 24616fdc2a6dSBarry Smith if (A->spd_set && A->spd) mumps->sym = 1; 24626fdc2a6dSBarry Smith else mumps->sym = 2; 246359ac8732SStefano Zampini #endif 2464450b117fSShri Abhyankar } 24652877fffaSHong Zhang 24662877fffaSHong Zhang mumps->isAIJ = PETSC_TRUE; 2467bf0cc555SLisandro Dalcin mumps->Destroy = B->ops->destroy; 24682877fffaSHong Zhang B->ops->destroy = MatDestroy_MUMPS; 24692877fffaSHong Zhang B->spptr = (void*)mumps; 24702205254eSKarl Rupp 2471f697e70eSHong Zhang ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 2472746480a1SHong Zhang 24732877fffaSHong Zhang *F = B; 24742877fffaSHong Zhang PetscFunctionReturn(0); 24752877fffaSHong Zhang } 24762877fffaSHong Zhang 2477bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI SBAIJ matrices */ 24782877fffaSHong Zhang #undef __FUNCT__ 2479bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_sbaij_mumps" 24808cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F) 24812877fffaSHong Zhang { 24822877fffaSHong Zhang Mat B; 24832877fffaSHong Zhang PetscErrorCode ierr; 24842877fffaSHong Zhang Mat_MUMPS *mumps; 2485ace3abfcSBarry Smith PetscBool isSeqSBAIJ; 24862877fffaSHong Zhang 24872877fffaSHong Zhang PetscFunctionBegin; 2488ce94432eSBarry Smith if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix"); 2489ce94432eSBarry 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"); 2490251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr); 24912877fffaSHong Zhang /* Create the factorization matrix */ 2492ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 24932877fffaSHong Zhang ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 24942877fffaSHong Zhang ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 2495b00a9115SJed Brown ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr); 2496bccb9932SShri Abhyankar if (isSeqSBAIJ) { 24970298fd71SBarry Smith ierr = MatSeqSBAIJSetPreallocation(B,1,0,NULL);CHKERRQ(ierr); 24982205254eSKarl Rupp 249916ebf90aSShri Abhyankar mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij; 2500dcd589f8SShri Abhyankar } else { 25010298fd71SBarry Smith ierr = MatMPISBAIJSetPreallocation(B,1,0,NULL,0,NULL);CHKERRQ(ierr); 25022205254eSKarl Rupp 2503bccb9932SShri Abhyankar mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij; 2504bccb9932SShri Abhyankar } 2505bccb9932SShri Abhyankar 250667877ebaSShri Abhyankar B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS; 2507bccb9932SShri Abhyankar B->ops->view = MatView_MUMPS; 250820be8e61SHong Zhang B->ops->getdiagonal = MatGetDiagonal_MUMPS; 25092205254eSKarl Rupp 2510bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 25115a05ddb0SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr); 25125a05ddb0SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorInvertSchurComplement_C",MatFactorInvertSchurComplement_MUMPS);CHKERRQ(ierr); 25135a05ddb0SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr); 25145a05ddb0SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSchurComplement_C",MatFactorGetSchurComplement_MUMPS);CHKERRQ(ierr); 25155a05ddb0SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplement_C",MatFactorSolveSchurComplement_MUMPS);CHKERRQ(ierr); 25165a05ddb0SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplementTranspose_C",MatFactorSolveSchurComplementTranspose_MUMPS);CHKERRQ(ierr); 2517b13644aeSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr); 2518b13644aeSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr); 2519b13644aeSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr); 2520b13644aeSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr); 2521ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr); 2522ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr); 2523ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr); 2524ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr); 25252205254eSKarl Rupp 2526f4762488SHong Zhang B->factortype = MAT_FACTOR_CHOLESKY; 252759ac8732SStefano Zampini #if defined(PETSC_USE_COMPLEX) 252859ac8732SStefano Zampini mumps->sym = 2; 252959ac8732SStefano Zampini #else 25306fdc2a6dSBarry Smith if (A->spd_set && A->spd) mumps->sym = 1; 25316fdc2a6dSBarry Smith else mumps->sym = 2; 253259ac8732SStefano Zampini #endif 2533a214ac2aSShri Abhyankar 2534bccb9932SShri Abhyankar mumps->isAIJ = PETSC_FALSE; 2535bf0cc555SLisandro Dalcin mumps->Destroy = B->ops->destroy; 2536f3c0ef26SHong Zhang B->ops->destroy = MatDestroy_MUMPS; 25372877fffaSHong Zhang B->spptr = (void*)mumps; 25382205254eSKarl Rupp 2539f697e70eSHong Zhang ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 2540746480a1SHong Zhang 25412877fffaSHong Zhang *F = B; 25422877fffaSHong Zhang PetscFunctionReturn(0); 25432877fffaSHong Zhang } 254497969023SHong Zhang 2545450b117fSShri Abhyankar #undef __FUNCT__ 2546bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_baij_mumps" 25478cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F) 254867877ebaSShri Abhyankar { 254967877ebaSShri Abhyankar Mat B; 255067877ebaSShri Abhyankar PetscErrorCode ierr; 255167877ebaSShri Abhyankar Mat_MUMPS *mumps; 2552ace3abfcSBarry Smith PetscBool isSeqBAIJ; 255367877ebaSShri Abhyankar 255467877ebaSShri Abhyankar PetscFunctionBegin; 255567877ebaSShri Abhyankar /* Create the factorization matrix */ 2556251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr); 2557ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 255867877ebaSShri Abhyankar ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 255967877ebaSShri Abhyankar ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 2560bccb9932SShri Abhyankar if (isSeqBAIJ) { 25610298fd71SBarry Smith ierr = MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,NULL);CHKERRQ(ierr); 2562bccb9932SShri Abhyankar } else { 25630298fd71SBarry Smith ierr = MatMPIBAIJSetPreallocation(B,A->rmap->bs,0,NULL,0,NULL);CHKERRQ(ierr); 2564bccb9932SShri Abhyankar } 2565450b117fSShri Abhyankar 2566b00a9115SJed Brown ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr); 2567450b117fSShri Abhyankar if (ftype == MAT_FACTOR_LU) { 2568450b117fSShri Abhyankar B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS; 2569450b117fSShri Abhyankar B->factortype = MAT_FACTOR_LU; 2570bccb9932SShri Abhyankar if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij; 2571bccb9932SShri Abhyankar else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij; 2572746480a1SHong Zhang mumps->sym = 0; 2573f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n"); 2574bccb9932SShri Abhyankar 2575450b117fSShri Abhyankar B->ops->view = MatView_MUMPS; 257620be8e61SHong Zhang B->ops->getdiagonal = MatGetDiagonal_MUMPS; 25772205254eSKarl Rupp 2578bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 25795a05ddb0SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr); 25805a05ddb0SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorInvertSchurComplement_C",MatFactorInvertSchurComplement_MUMPS);CHKERRQ(ierr); 25815a05ddb0SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr); 25825a05ddb0SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSchurComplement_C",MatFactorGetSchurComplement_MUMPS);CHKERRQ(ierr); 25835a05ddb0SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplement_C",MatFactorSolveSchurComplement_MUMPS);CHKERRQ(ierr); 25845a05ddb0SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplementTranspose_C",MatFactorSolveSchurComplementTranspose_MUMPS);CHKERRQ(ierr); 2585bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr); 2586bc6112feSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr); 2587bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr); 2588bc6112feSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr); 2589ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr); 2590ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr); 2591ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr); 2592ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr); 2593450b117fSShri Abhyankar 2594450b117fSShri Abhyankar mumps->isAIJ = PETSC_TRUE; 2595bf0cc555SLisandro Dalcin mumps->Destroy = B->ops->destroy; 2596450b117fSShri Abhyankar B->ops->destroy = MatDestroy_MUMPS; 2597450b117fSShri Abhyankar B->spptr = (void*)mumps; 25982205254eSKarl Rupp 2599f697e70eSHong Zhang ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 2600746480a1SHong Zhang 2601450b117fSShri Abhyankar *F = B; 2602450b117fSShri Abhyankar PetscFunctionReturn(0); 2603450b117fSShri Abhyankar } 260442c9c57cSBarry Smith 260542c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mumps(Mat,MatFactorType,Mat*); 260642c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_baij_mumps(Mat,MatFactorType,Mat*); 260742c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat,MatFactorType,Mat*); 260842c9c57cSBarry Smith 260942c9c57cSBarry Smith #undef __FUNCT__ 261042c9c57cSBarry Smith #define __FUNCT__ "MatSolverPackageRegister_MUMPS" 261129b38603SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_MUMPS(void) 261242c9c57cSBarry Smith { 261342c9c57cSBarry Smith PetscErrorCode ierr; 261442c9c57cSBarry Smith 261542c9c57cSBarry Smith PetscFunctionBegin; 261642c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr); 261742c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr); 261842c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr); 261942c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr); 262042c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPISBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr); 262142c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr); 262242c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr); 262342c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr); 262442c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr); 262542c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr); 262642c9c57cSBarry Smith PetscFunctionReturn(0); 262742c9c57cSBarry Smith } 262842c9c57cSBarry Smith 2629