148f88336SPieter Ghysels #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 27d6ea485SPieter Ghysels #include <../src/mat/impls/aij/mpi/mpiaij.h> 37d6ea485SPieter Ghysels #include <StrumpackSparseSolver.h> 47d6ea485SPieter Ghysels 59371c9d4SSatish Balay static PetscErrorCode MatGetDiagonal_STRUMPACK(Mat A, Vec v) { 67d6ea485SPieter Ghysels PetscFunctionBegin; 77d6ea485SPieter Ghysels SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Mat type: STRUMPACK factor"); 87d6ea485SPieter Ghysels PetscFunctionReturn(0); 97d6ea485SPieter Ghysels } 107d6ea485SPieter Ghysels 119371c9d4SSatish Balay static PetscErrorCode MatDestroy_STRUMPACK(Mat A) { 120d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)A->spptr; 137d6ea485SPieter Ghysels PetscBool flg; 147d6ea485SPieter Ghysels 157d6ea485SPieter Ghysels PetscFunctionBegin; 167d6ea485SPieter Ghysels /* Deallocate STRUMPACK storage */ 17e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_destroy", STRUMPACK_destroy(S)); 189566063dSJacob Faibussowitsch PetscCall(PetscFree(A->spptr)); 199566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &flg)); 207d6ea485SPieter Ghysels if (flg) { 219566063dSJacob Faibussowitsch PetscCall(MatDestroy_SeqAIJ(A)); 227d6ea485SPieter Ghysels } else { 239566063dSJacob Faibussowitsch PetscCall(MatDestroy_MPIAIJ(A)); 247d6ea485SPieter Ghysels } 25575704cbSPieter Ghysels 26575704cbSPieter Ghysels /* clear composed functions */ 279566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorGetSolverType_C", NULL)); 289566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetReordering_C", NULL)); 299566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetColPerm_C", NULL)); 309566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetHSSRelTol_C", NULL)); 319566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetHSSAbsTol_C", NULL)); 329566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetHSSMaxRank_C", NULL)); 339566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetHSSLeafSize_C", NULL)); 349566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetHSSMinSepSize_C", NULL)); 35575704cbSPieter Ghysels 36575704cbSPieter Ghysels PetscFunctionReturn(0); 37575704cbSPieter Ghysels } 38575704cbSPieter Ghysels 399371c9d4SSatish Balay static PetscErrorCode MatSTRUMPACKSetReordering_STRUMPACK(Mat F, MatSTRUMPACKReordering reordering) { 400d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->spptr; 41575704cbSPieter Ghysels 42575704cbSPieter Ghysels PetscFunctionBegin; 43e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_reordering_method", STRUMPACK_set_reordering_method(*S, (STRUMPACK_REORDERING_STRATEGY)reordering)); 4434f43fa5SPieter Ghysels PetscFunctionReturn(0); 4534f43fa5SPieter Ghysels } 46e5a36eccSStefano Zampini 47542aee0fSPieter Ghysels /*@ 4834f43fa5SPieter Ghysels MatSTRUMPACKSetReordering - Set STRUMPACK fill-reducing reordering 4934f43fa5SPieter Ghysels 5034f43fa5SPieter Ghysels Input Parameters: 5111a5261eSBarry Smith + F - the factored matrix obtained by calling `MatGetFactor(`) from PETSc-STRUMPACK interface 5234f43fa5SPieter Ghysels - reordering - the code to be used to find the fill-reducing reordering 5334f43fa5SPieter Ghysels Possible values: NATURAL=0 METIS=1 PARMETIS=2 SCOTCH=3 PTSCOTCH=4 RCM=5 5434f43fa5SPieter Ghysels 5511a5261eSBarry Smith Options Database Key: 5634f43fa5SPieter Ghysels . -mat_strumpack_reordering <METIS> - Sparsity reducing matrix reordering (choose one of) NATURAL METIS PARMETIS SCOTCH PTSCOTCH RCM (None) 5734f43fa5SPieter Ghysels 5834f43fa5SPieter Ghysels Level: beginner 5934f43fa5SPieter Ghysels 6034f43fa5SPieter Ghysels References: 61606c0280SSatish Balay . * - STRUMPACK manual 6234f43fa5SPieter Ghysels 63db781477SPatrick Sanan .seealso: `MatGetFactor()` 64542aee0fSPieter Ghysels @*/ 659371c9d4SSatish Balay PetscErrorCode MatSTRUMPACKSetReordering(Mat F, MatSTRUMPACKReordering reordering) { 6634f43fa5SPieter Ghysels PetscFunctionBegin; 6734f43fa5SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 6834f43fa5SPieter Ghysels PetscValidLogicalCollectiveEnum(F, reordering, 2); 69cac4c232SBarry Smith PetscTryMethod(F, "MatSTRUMPACKSetReordering_C", (Mat, MatSTRUMPACKReordering), (F, reordering)); 70575704cbSPieter Ghysels PetscFunctionReturn(0); 71575704cbSPieter Ghysels } 72575704cbSPieter Ghysels 739371c9d4SSatish Balay static PetscErrorCode MatSTRUMPACKSetColPerm_STRUMPACK(Mat F, PetscBool cperm) { 7434f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->spptr; 7534f43fa5SPieter Ghysels 7634f43fa5SPieter Ghysels PetscFunctionBegin; 77e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_set_mc64job", STRUMPACK_set_mc64job(*S, cperm ? 5 : 0)); 7834f43fa5SPieter Ghysels PetscFunctionReturn(0); 7934f43fa5SPieter Ghysels } 80e5a36eccSStefano Zampini 81575704cbSPieter Ghysels /*@ 8234f43fa5SPieter Ghysels MatSTRUMPACKSetColPerm - Set whether STRUMPACK should try to permute the columns of the matrix in order to get a nonzero diagonal 83147403d9SBarry Smith 8411a5261eSBarry Smith Logically Collective on F 85575704cbSPieter Ghysels 86575704cbSPieter Ghysels Input Parameters: 8711a5261eSBarry Smith + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 8811a5261eSBarry Smith - cperm - `PETSC_TRUE` to permute (internally) the columns of the matrix 89575704cbSPieter Ghysels 9011a5261eSBarry Smith Options Database Key: 91147403d9SBarry Smith . -mat_strumpack_colperm <cperm> - true to use the permutation 92575704cbSPieter Ghysels 93575704cbSPieter Ghysels Level: beginner 94575704cbSPieter Ghysels 95575704cbSPieter Ghysels References: 96606c0280SSatish Balay . * - STRUMPACK manual 97575704cbSPieter Ghysels 98db781477SPatrick Sanan .seealso: `MatGetFactor()` 99575704cbSPieter Ghysels @*/ 1009371c9d4SSatish Balay PetscErrorCode MatSTRUMPACKSetColPerm(Mat F, PetscBool cperm) { 101575704cbSPieter Ghysels PetscFunctionBegin; 102575704cbSPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 10334f43fa5SPieter Ghysels PetscValidLogicalCollectiveBool(F, cperm, 2); 104cac4c232SBarry Smith PetscTryMethod(F, "MatSTRUMPACKSetColPerm_C", (Mat, PetscBool), (F, cperm)); 105575704cbSPieter Ghysels PetscFunctionReturn(0); 106575704cbSPieter Ghysels } 107575704cbSPieter Ghysels 1089371c9d4SSatish Balay static PetscErrorCode MatSTRUMPACKSetHSSRelTol_STRUMPACK(Mat F, PetscReal rtol) { 10934f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->spptr; 11034f43fa5SPieter Ghysels 11134f43fa5SPieter Ghysels PetscFunctionBegin; 112e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_set_HSS_rel_tol", STRUMPACK_set_HSS_rel_tol(*S, rtol)); 11334f43fa5SPieter Ghysels PetscFunctionReturn(0); 11434f43fa5SPieter Ghysels } 115e5a36eccSStefano Zampini 11634f43fa5SPieter Ghysels /*@ 11734f43fa5SPieter Ghysels MatSTRUMPACKSetHSSRelTol - Set STRUMPACK relative tolerance for HSS compression 118147403d9SBarry Smith 11911a5261eSBarry Smith Logically Collective on F 12034f43fa5SPieter Ghysels 12134f43fa5SPieter Ghysels Input Parameters: 12211a5261eSBarry Smith + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 12334f43fa5SPieter Ghysels - rtol - relative compression tolerance 12434f43fa5SPieter Ghysels 12511a5261eSBarry Smith Options Database Key: 12634f43fa5SPieter Ghysels . -mat_strumpack_hss_rel_tol <1e-2> - Relative compression tolerance (None) 12734f43fa5SPieter Ghysels 12834f43fa5SPieter Ghysels Level: beginner 12934f43fa5SPieter Ghysels 13034f43fa5SPieter Ghysels References: 131606c0280SSatish Balay . * - STRUMPACK manual 13234f43fa5SPieter Ghysels 133db781477SPatrick Sanan .seealso: `MatGetFactor()` 13434f43fa5SPieter Ghysels @*/ 1359371c9d4SSatish Balay PetscErrorCode MatSTRUMPACKSetHSSRelTol(Mat F, PetscReal rtol) { 13634f43fa5SPieter Ghysels PetscFunctionBegin; 13734f43fa5SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 13834f43fa5SPieter Ghysels PetscValidLogicalCollectiveReal(F, rtol, 2); 139cac4c232SBarry Smith PetscTryMethod(F, "MatSTRUMPACKSetHSSRelTol_C", (Mat, PetscReal), (F, rtol)); 14034f43fa5SPieter Ghysels PetscFunctionReturn(0); 14134f43fa5SPieter Ghysels } 14234f43fa5SPieter Ghysels 1439371c9d4SSatish Balay static PetscErrorCode MatSTRUMPACKSetHSSAbsTol_STRUMPACK(Mat F, PetscReal atol) { 14434f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->spptr; 14534f43fa5SPieter Ghysels 14634f43fa5SPieter Ghysels PetscFunctionBegin; 147e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_set_HSS_abs_tol", STRUMPACK_set_HSS_abs_tol(*S, atol)); 14834f43fa5SPieter Ghysels PetscFunctionReturn(0); 14934f43fa5SPieter Ghysels } 150e5a36eccSStefano Zampini 15134f43fa5SPieter Ghysels /*@ 15234f43fa5SPieter Ghysels MatSTRUMPACKSetHSSAbsTol - Set STRUMPACK absolute tolerance for HSS compression 153147403d9SBarry Smith 15411a5261eSBarry Smith Logically Collective on F 15534f43fa5SPieter Ghysels 15634f43fa5SPieter Ghysels Input Parameters: 15711a5261eSBarry Smith + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 15834f43fa5SPieter Ghysels - atol - absolute compression tolerance 15934f43fa5SPieter Ghysels 16011a5261eSBarry Smith Options Database Key: 16134f43fa5SPieter Ghysels . -mat_strumpack_hss_abs_tol <1e-8> - Absolute compression tolerance (None) 16234f43fa5SPieter Ghysels 16334f43fa5SPieter Ghysels Level: beginner 16434f43fa5SPieter Ghysels 16534f43fa5SPieter Ghysels References: 166606c0280SSatish Balay . * - STRUMPACK manual 16734f43fa5SPieter Ghysels 168db781477SPatrick Sanan .seealso: `MatGetFactor()` 16934f43fa5SPieter Ghysels @*/ 1709371c9d4SSatish Balay PetscErrorCode MatSTRUMPACKSetHSSAbsTol(Mat F, PetscReal atol) { 17134f43fa5SPieter Ghysels PetscFunctionBegin; 17234f43fa5SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 17334f43fa5SPieter Ghysels PetscValidLogicalCollectiveReal(F, atol, 2); 174cac4c232SBarry Smith PetscTryMethod(F, "MatSTRUMPACKSetHSSAbsTol_C", (Mat, PetscReal), (F, atol)); 17534f43fa5SPieter Ghysels PetscFunctionReturn(0); 17634f43fa5SPieter Ghysels } 17734f43fa5SPieter Ghysels 1789371c9d4SSatish Balay static PetscErrorCode MatSTRUMPACKSetHSSMaxRank_STRUMPACK(Mat F, PetscInt hssmaxrank) { 17934f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->spptr; 18034f43fa5SPieter Ghysels 18134f43fa5SPieter Ghysels PetscFunctionBegin; 182e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_set_HSS_max_rank", STRUMPACK_set_HSS_max_rank(*S, hssmaxrank)); 18334f43fa5SPieter Ghysels PetscFunctionReturn(0); 18434f43fa5SPieter Ghysels } 185e5a36eccSStefano Zampini 18634f43fa5SPieter Ghysels /*@ 18734f43fa5SPieter Ghysels MatSTRUMPACKSetHSSMaxRank - Set STRUMPACK maximum HSS rank 188147403d9SBarry Smith 18911a5261eSBarry Smith Logically Collective on F 19034f43fa5SPieter Ghysels 19134f43fa5SPieter Ghysels Input Parameters: 19211a5261eSBarry Smith + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 19334f43fa5SPieter Ghysels - hssmaxrank - maximum rank used in low-rank approximation 19434f43fa5SPieter Ghysels 19511a5261eSBarry Smith Options Database Key: 19634f43fa5SPieter Ghysels . -mat_strumpack_max_rank - Maximum rank in HSS compression, when using pctype ilu (None) 19734f43fa5SPieter Ghysels 19834f43fa5SPieter Ghysels Level: beginner 19934f43fa5SPieter Ghysels 20034f43fa5SPieter Ghysels References: 201606c0280SSatish Balay . * - STRUMPACK manual 20234f43fa5SPieter Ghysels 203db781477SPatrick Sanan .seealso: `MatGetFactor()` 20434f43fa5SPieter Ghysels @*/ 2059371c9d4SSatish Balay PetscErrorCode MatSTRUMPACKSetHSSMaxRank(Mat F, PetscInt hssmaxrank) { 20634f43fa5SPieter Ghysels PetscFunctionBegin; 20734f43fa5SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 20834f43fa5SPieter Ghysels PetscValidLogicalCollectiveInt(F, hssmaxrank, 2); 209cac4c232SBarry Smith PetscTryMethod(F, "MatSTRUMPACKSetHSSMaxRank_C", (Mat, PetscInt), (F, hssmaxrank)); 21034f43fa5SPieter Ghysels PetscFunctionReturn(0); 21134f43fa5SPieter Ghysels } 21234f43fa5SPieter Ghysels 2139371c9d4SSatish Balay static PetscErrorCode MatSTRUMPACKSetHSSLeafSize_STRUMPACK(Mat F, PetscInt leaf_size) { 214a36bf211SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->spptr; 215a36bf211SPieter Ghysels 216a36bf211SPieter Ghysels PetscFunctionBegin; 217e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_set_HSS_leaf_size", STRUMPACK_set_HSS_leaf_size(*S, leaf_size)); 218a36bf211SPieter Ghysels PetscFunctionReturn(0); 219a36bf211SPieter Ghysels } 220e5a36eccSStefano Zampini 221a36bf211SPieter Ghysels /*@ 222a36bf211SPieter Ghysels MatSTRUMPACKSetHSSLeafSize - Set STRUMPACK HSS leaf size 223147403d9SBarry Smith 22411a5261eSBarry Smith Logically Collective on F 225a36bf211SPieter Ghysels 226a36bf211SPieter Ghysels Input Parameters: 22711a5261eSBarry Smith + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 228a36bf211SPieter Ghysels - leaf_size - Size of diagonal blocks in HSS approximation 229a36bf211SPieter Ghysels 23011a5261eSBarry Smith Options Database Key: 231a36bf211SPieter Ghysels . -mat_strumpack_leaf_size - Size of diagonal blocks in HSS approximation, when using pctype ilu (None) 232a36bf211SPieter Ghysels 233a36bf211SPieter Ghysels Level: beginner 234a36bf211SPieter Ghysels 235a36bf211SPieter Ghysels References: 236606c0280SSatish Balay . * - STRUMPACK manual 237a36bf211SPieter Ghysels 238db781477SPatrick Sanan .seealso: `MatGetFactor()` 239a36bf211SPieter Ghysels @*/ 2409371c9d4SSatish Balay PetscErrorCode MatSTRUMPACKSetHSSLeafSize(Mat F, PetscInt leaf_size) { 241a36bf211SPieter Ghysels PetscFunctionBegin; 242a36bf211SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 243a36bf211SPieter Ghysels PetscValidLogicalCollectiveInt(F, leaf_size, 2); 244cac4c232SBarry Smith PetscTryMethod(F, "MatSTRUMPACKSetHSSLeafSize_C", (Mat, PetscInt), (F, leaf_size)); 245a36bf211SPieter Ghysels PetscFunctionReturn(0); 246a36bf211SPieter Ghysels } 247a36bf211SPieter Ghysels 2489371c9d4SSatish Balay static PetscErrorCode MatSTRUMPACKSetHSSMinSepSize_STRUMPACK(Mat F, PetscInt hssminsize) { 24934f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->spptr; 25034f43fa5SPieter Ghysels 25134f43fa5SPieter Ghysels PetscFunctionBegin; 252e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_set_HSS_min_sep_size", STRUMPACK_set_HSS_min_sep_size(*S, hssminsize)); 25334f43fa5SPieter Ghysels PetscFunctionReturn(0); 25434f43fa5SPieter Ghysels } 255e5a36eccSStefano Zampini 256291d0ed5SPieter Ghysels /*@ 257291d0ed5SPieter Ghysels MatSTRUMPACKSetHSSMinSepSize - Set STRUMPACK minimum separator size for low-rank approximation 258147403d9SBarry Smith 25911a5261eSBarry Smith Logically Collective on F 260291d0ed5SPieter Ghysels 261291d0ed5SPieter Ghysels Input Parameters: 26211a5261eSBarry Smith + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 263291d0ed5SPieter Ghysels - hssminsize - minimum dense matrix size for low-rank approximation 264291d0ed5SPieter Ghysels 26511a5261eSBarry Smith Options Database Key: 266147403d9SBarry Smith . -mat_strumpack_hss_min_sep_size <hssminsize> - set the minimum separator size 267291d0ed5SPieter Ghysels 268291d0ed5SPieter Ghysels Level: beginner 269291d0ed5SPieter Ghysels 270291d0ed5SPieter Ghysels References: 271606c0280SSatish Balay . * - STRUMPACK manual 272291d0ed5SPieter Ghysels 273db781477SPatrick Sanan .seealso: `MatGetFactor()` 274291d0ed5SPieter Ghysels @*/ 2759371c9d4SSatish Balay PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat F, PetscInt hssminsize) { 276291d0ed5SPieter Ghysels PetscFunctionBegin; 277291d0ed5SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 278291d0ed5SPieter Ghysels PetscValidLogicalCollectiveInt(F, hssminsize, 2); 279cac4c232SBarry Smith PetscTryMethod(F, "MatSTRUMPACKSetHSSMinSepSize_C", (Mat, PetscInt), (F, hssminsize)); 280291d0ed5SPieter Ghysels PetscFunctionReturn(0); 281291d0ed5SPieter Ghysels } 282291d0ed5SPieter Ghysels 2839371c9d4SSatish Balay static PetscErrorCode MatSolve_STRUMPACK(Mat A, Vec b_mpi, Vec x) { 2840d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)A->spptr; 2857d6ea485SPieter Ghysels STRUMPACK_RETURN_CODE sp_err; 2867d6ea485SPieter Ghysels const PetscScalar *bptr; 2877d6ea485SPieter Ghysels PetscScalar *xptr; 2887d6ea485SPieter Ghysels 2897d6ea485SPieter Ghysels PetscFunctionBegin; 2909566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &xptr)); 2919566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(b_mpi, &bptr)); 2920d08a34dSPieter Ghysels 293e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_solve", sp_err = STRUMPACK_solve(*S, (PetscScalar *)bptr, xptr, 0)); 2940d08a34dSPieter Ghysels switch (sp_err) { 2950d08a34dSPieter Ghysels case STRUMPACK_SUCCESS: break; 2969371c9d4SSatish Balay case STRUMPACK_MATRIX_NOT_SET: { 2979371c9d4SSatish Balay SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix was not set"); 2989371c9d4SSatish Balay break; 2999371c9d4SSatish Balay } 3009371c9d4SSatish Balay case STRUMPACK_REORDERING_ERROR: { 3019371c9d4SSatish Balay SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix reordering failed"); 3029371c9d4SSatish Balay break; 3039371c9d4SSatish Balay } 3040d08a34dSPieter Ghysels default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: solve failed"); 3057d6ea485SPieter Ghysels } 3069566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &xptr)); 3079566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(b_mpi, &bptr)); 3087d6ea485SPieter Ghysels PetscFunctionReturn(0); 3097d6ea485SPieter Ghysels } 3107d6ea485SPieter Ghysels 3119371c9d4SSatish Balay static PetscErrorCode MatMatSolve_STRUMPACK(Mat A, Mat B_mpi, Mat X) { 3127d6ea485SPieter Ghysels PetscBool flg; 3137d6ea485SPieter Ghysels 3147d6ea485SPieter Ghysels PetscFunctionBegin; 3159566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)B_mpi, &flg, MATSEQDENSE, MATMPIDENSE, NULL)); 3165f80ce2aSJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Matrix B must be MATDENSE matrix"); 3179566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)X, &flg, MATSEQDENSE, MATMPIDENSE, NULL)); 3185f80ce2aSJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Matrix X must be MATDENSE matrix"); 3197d6ea485SPieter Ghysels SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "MatMatSolve_STRUMPACK() is not implemented yet"); 3207d6ea485SPieter Ghysels PetscFunctionReturn(0); 3217d6ea485SPieter Ghysels } 3227d6ea485SPieter Ghysels 3239371c9d4SSatish Balay static PetscErrorCode MatView_Info_STRUMPACK(Mat A, PetscViewer viewer) { 324ad0c5e61SPieter Ghysels PetscFunctionBegin; 325ad0c5e61SPieter Ghysels /* check if matrix is strumpack type */ 326ad0c5e61SPieter Ghysels if (A->ops->solve != MatSolve_STRUMPACK) PetscFunctionReturn(0); 3279566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "STRUMPACK sparse solver!\n")); 328ad0c5e61SPieter Ghysels PetscFunctionReturn(0); 329ad0c5e61SPieter Ghysels } 330ad0c5e61SPieter Ghysels 3319371c9d4SSatish Balay static PetscErrorCode MatView_STRUMPACK(Mat A, PetscViewer viewer) { 332ad0c5e61SPieter Ghysels PetscBool iascii; 333ad0c5e61SPieter Ghysels PetscViewerFormat format; 334ad0c5e61SPieter Ghysels 335ad0c5e61SPieter Ghysels PetscFunctionBegin; 3369566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 337ad0c5e61SPieter Ghysels if (iascii) { 3389566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 3399566063dSJacob Faibussowitsch if (format == PETSC_VIEWER_ASCII_INFO) PetscCall(MatView_Info_STRUMPACK(A, viewer)); 340ad0c5e61SPieter Ghysels } 341ad0c5e61SPieter Ghysels PetscFunctionReturn(0); 342ad0c5e61SPieter Ghysels } 3437d6ea485SPieter Ghysels 3449371c9d4SSatish Balay static PetscErrorCode MatLUFactorNumeric_STRUMPACK(Mat F, Mat A, const MatFactorInfo *info) { 3450d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->spptr; 3467d6ea485SPieter Ghysels STRUMPACK_RETURN_CODE sp_err; 3477d6ea485SPieter Ghysels Mat_SeqAIJ *A_d, *A_o; 3487d6ea485SPieter Ghysels Mat_MPIAIJ *mat; 3490d08a34dSPieter Ghysels PetscInt M = A->rmap->N, m = A->rmap->n; 3507d6ea485SPieter Ghysels PetscBool flg; 3517d6ea485SPieter Ghysels 3527d6ea485SPieter Ghysels PetscFunctionBegin; 3539566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIAIJ, &flg)); 3547d6ea485SPieter Ghysels if (flg) { /* A is MATMPIAIJ */ 3557d6ea485SPieter Ghysels mat = (Mat_MPIAIJ *)A->data; 3567d6ea485SPieter Ghysels A_d = (Mat_SeqAIJ *)(mat->A)->data; 3577d6ea485SPieter Ghysels A_o = (Mat_SeqAIJ *)(mat->B)->data; 358e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_set_MPIAIJ_matrix", STRUMPACK_set_MPIAIJ_matrix(*S, &m, A_d->i, A_d->j, A_d->a, A_o->i, A_o->j, A_o->a, mat->garray)); 3597d6ea485SPieter Ghysels } else { /* A is MATSEQAIJ */ 3607d6ea485SPieter Ghysels A_d = (Mat_SeqAIJ *)A->data; 361e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_set_csr_matrix", STRUMPACK_set_csr_matrix(*S, &M, A_d->i, A_d->j, A_d->a, 0)); 3627d6ea485SPieter Ghysels } 3637d6ea485SPieter Ghysels 3647d6ea485SPieter Ghysels /* Reorder and Factor the matrix. */ 3657d6ea485SPieter Ghysels /* TODO figure out how to avoid reorder if the matrix values changed, but the pattern remains the same. */ 366e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_reorder", sp_err = STRUMPACK_reorder(*S)); 367e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_factor", sp_err = STRUMPACK_factor(*S)); 3680d08a34dSPieter Ghysels switch (sp_err) { 3690d08a34dSPieter Ghysels case STRUMPACK_SUCCESS: break; 3709371c9d4SSatish Balay case STRUMPACK_MATRIX_NOT_SET: { 3719371c9d4SSatish Balay SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix was not set"); 3729371c9d4SSatish Balay break; 3739371c9d4SSatish Balay } 3749371c9d4SSatish Balay case STRUMPACK_REORDERING_ERROR: { 3759371c9d4SSatish Balay SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix reordering failed"); 3769371c9d4SSatish Balay break; 3779371c9d4SSatish Balay } 3780d08a34dSPieter Ghysels default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: factorization failed"); 3797d6ea485SPieter Ghysels } 380cb250fa3SPieter Ghysels F->assembled = PETSC_TRUE; 381cb250fa3SPieter Ghysels F->preallocated = PETSC_TRUE; 3827d6ea485SPieter Ghysels PetscFunctionReturn(0); 3837d6ea485SPieter Ghysels } 3847d6ea485SPieter Ghysels 3859371c9d4SSatish Balay static PetscErrorCode MatLUFactorSymbolic_STRUMPACK(Mat F, Mat A, IS r, IS c, const MatFactorInfo *info) { 38626cc229bSBarry Smith STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->spptr; 38726cc229bSBarry Smith PetscBool flg, set; 38826cc229bSBarry Smith PetscReal ctol; 38926cc229bSBarry Smith PetscInt hssminsize, max_rank, leaf_size; 39026cc229bSBarry Smith STRUMPACK_REORDERING_STRATEGY ndcurrent, ndvalue; 39126cc229bSBarry Smith STRUMPACK_KRYLOV_SOLVER itcurrent, itsolver; 39226cc229bSBarry Smith const char *const STRUMPACKNDTypes[] = {"NATURAL", "METIS", "PARMETIS", "SCOTCH", "PTSCOTCH", "RCM", "STRUMPACKNDTypes", "", 0}; 39326cc229bSBarry Smith const char *const SolverTypes[] = {"AUTO", "NONE", "REFINE", "PREC_GMRES", "GMRES", "PREC_BICGSTAB", "BICGSTAB", "SolverTypes", "", 0}; 39426cc229bSBarry Smith 3957d6ea485SPieter Ghysels PetscFunctionBegin; 39626cc229bSBarry Smith /* Set options to F */ 39726cc229bSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)F), ((PetscObject)F)->prefix, "STRUMPACK Options", "Mat"); 398e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_HSS_rel_tol", ctol = (PetscReal)STRUMPACK_HSS_rel_tol(*S)); 39926cc229bSBarry Smith PetscCall(PetscOptionsReal("-mat_strumpack_hss_rel_tol", "Relative compression tolerance", "None", ctol, &ctol, &set)); 400e77caa6dSBarry Smith if (set) PetscStackCallExternalVoid("STRUMPACK_set_HSS_rel_tol", STRUMPACK_set_HSS_rel_tol(*S, (double)ctol)); 40126cc229bSBarry Smith 402e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_HSS_abs_tol", ctol = (PetscReal)STRUMPACK_HSS_abs_tol(*S)); 40326cc229bSBarry Smith PetscCall(PetscOptionsReal("-mat_strumpack_hss_abs_tol", "Absolute compression tolerance", "None", ctol, &ctol, &set)); 404e77caa6dSBarry Smith if (set) PetscStackCallExternalVoid("STRUMPACK_set_HSS_abs_tol", STRUMPACK_set_HSS_abs_tol(*S, (double)ctol)); 40526cc229bSBarry Smith 406e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_mc64job", flg = (STRUMPACK_mc64job(*S) == 0) ? PETSC_FALSE : PETSC_TRUE); 40726cc229bSBarry Smith PetscCall(PetscOptionsBool("-mat_strumpack_colperm", "Find a col perm to get nonzero diagonal", "None", flg, &flg, &set)); 408e77caa6dSBarry Smith if (set) PetscStackCallExternalVoid("STRUMPACK_set_mc64job", STRUMPACK_set_mc64job(*S, flg ? 5 : 0)); 40926cc229bSBarry Smith 410e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_HSS_min_sep_size", hssminsize = (PetscInt)STRUMPACK_HSS_min_sep_size(*S)); 41126cc229bSBarry Smith PetscCall(PetscOptionsInt("-mat_strumpack_hss_min_sep_size", "Minimum size of separator for HSS compression", "None", hssminsize, &hssminsize, &set)); 412e77caa6dSBarry Smith if (set) PetscStackCallExternalVoid("STRUMPACK_set_HSS_min_sep_size", STRUMPACK_set_HSS_min_sep_size(*S, (int)hssminsize)); 41326cc229bSBarry Smith 414e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_HSS_max_rank", max_rank = (PetscInt)STRUMPACK_HSS_max_rank(*S)); 41526cc229bSBarry Smith PetscCall(PetscOptionsInt("-mat_strumpack_max_rank", "Maximum rank in HSS approximation", "None", max_rank, &max_rank, &set)); 416e77caa6dSBarry Smith if (set) PetscStackCallExternalVoid("STRUMPACK_set_HSS_max_rank", STRUMPACK_set_HSS_max_rank(*S, (int)max_rank)); 41726cc229bSBarry Smith 418e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_HSS_leaf_size", leaf_size = (PetscInt)STRUMPACK_HSS_leaf_size(*S)); 41926cc229bSBarry Smith PetscCall(PetscOptionsInt("-mat_strumpack_leaf_size", "Size of diagonal blocks in HSS approximation", "None", leaf_size, &leaf_size, &set)); 420e77caa6dSBarry Smith if (set) PetscStackCallExternalVoid("STRUMPACK_set_HSS_leaf_size", STRUMPACK_set_HSS_leaf_size(*S, (int)leaf_size)); 42126cc229bSBarry Smith 422e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_reordering_method", ndcurrent = STRUMPACK_reordering_method(*S)); 42326cc229bSBarry Smith PetscCall(PetscOptionsEnum("-mat_strumpack_reordering", "Sparsity reducing matrix reordering", "None", STRUMPACKNDTypes, (PetscEnum)ndcurrent, (PetscEnum *)&ndvalue, &set)); 424e77caa6dSBarry Smith if (set) PetscStackCallExternalVoid("STRUMPACK_set_reordering_method", STRUMPACK_set_reordering_method(*S, ndvalue)); 42526cc229bSBarry Smith 42626cc229bSBarry Smith /* Disable the outer iterative solver from STRUMPACK. */ 42726cc229bSBarry Smith /* When STRUMPACK is used as a direct solver, it will by default do iterative refinement. */ 42826cc229bSBarry Smith /* When STRUMPACK is used as an approximate factorization preconditioner (by enabling */ 42926cc229bSBarry Smith /* low-rank compression), it will use it's own preconditioned GMRES. Here we can disable */ 43026cc229bSBarry Smith /* the outer iterative solver, as PETSc uses STRUMPACK from within a KSP. */ 431e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_set_Krylov_solver", STRUMPACK_set_Krylov_solver(*S, STRUMPACK_DIRECT)); 43226cc229bSBarry Smith 433e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_Krylov_solver", itcurrent = STRUMPACK_Krylov_solver(*S)); 43426cc229bSBarry Smith PetscCall(PetscOptionsEnum("-mat_strumpack_iterative_solver", "Select iterative solver from STRUMPACK", "None", SolverTypes, (PetscEnum)itcurrent, (PetscEnum *)&itsolver, &set)); 435e77caa6dSBarry Smith if (set) PetscStackCallExternalVoid("STRUMPACK_set_Krylov_solver", STRUMPACK_set_Krylov_solver(*S, itsolver)); 43626cc229bSBarry Smith PetscOptionsEnd(); 43726cc229bSBarry Smith 4387d6ea485SPieter Ghysels F->ops->lufactornumeric = MatLUFactorNumeric_STRUMPACK; 4397d6ea485SPieter Ghysels F->ops->solve = MatSolve_STRUMPACK; 4407d6ea485SPieter Ghysels F->ops->matsolve = MatMatSolve_STRUMPACK; 4417d6ea485SPieter Ghysels PetscFunctionReturn(0); 4427d6ea485SPieter Ghysels } 4437d6ea485SPieter Ghysels 4449371c9d4SSatish Balay static PetscErrorCode MatFactorGetSolverType_aij_strumpack(Mat A, MatSolverType *type) { 4457d6ea485SPieter Ghysels PetscFunctionBegin; 4467d6ea485SPieter Ghysels *type = MATSOLVERSTRUMPACK; 4477d6ea485SPieter Ghysels PetscFunctionReturn(0); 4487d6ea485SPieter Ghysels } 4497d6ea485SPieter Ghysels 450575704cbSPieter Ghysels /*MC 451575704cbSPieter Ghysels MATSOLVERSSTRUMPACK = "strumpack" - A solver package providing a direct sparse solver (PCLU) 452575704cbSPieter Ghysels and a preconditioner (PCILU) using low-rank compression via the external package STRUMPACK. 453575704cbSPieter Ghysels 454575704cbSPieter Ghysels Consult the STRUMPACK-sparse manual for more info. 455575704cbSPieter Ghysels 456575704cbSPieter Ghysels Use 457575704cbSPieter Ghysels ./configure --download-strumpack 458575704cbSPieter Ghysels to have PETSc installed with STRUMPACK 459575704cbSPieter Ghysels 460575704cbSPieter Ghysels Use 4613ca39a21SBarry Smith -pc_type lu -pc_factor_mat_solver_type strumpack 462575704cbSPieter Ghysels to use this as an exact (direct) solver, use 4633ca39a21SBarry Smith -pc_type ilu -pc_factor_mat_solver_type strumpack 464575704cbSPieter Ghysels to enable low-rank compression (i.e, use as a preconditioner). 465575704cbSPieter Ghysels 466575704cbSPieter Ghysels Works with AIJ matrices 467575704cbSPieter Ghysels 468575704cbSPieter Ghysels Options Database Keys: 46967b8a455SSatish Balay + -mat_strumpack_verbose - verbose info 47034f43fa5SPieter Ghysels . -mat_strumpack_hss_rel_tol <1e-2> - Relative compression tolerance (None) 47134f43fa5SPieter Ghysels . -mat_strumpack_hss_abs_tol <1e-8> - Absolute compression tolerance (None) 472575704cbSPieter Ghysels . -mat_strumpack_colperm <TRUE> - Permute matrix to make diagonal nonzeros (None) 47334f43fa5SPieter Ghysels . -mat_strumpack_hss_min_sep_size <256> - Minimum size of separator for HSS compression (None) 47434f43fa5SPieter Ghysels . -mat_strumpack_max_rank - Maximum rank in HSS compression, when using pctype ilu (None) 475a36bf211SPieter Ghysels . -mat_strumpack_leaf_size - Size of diagonal blocks in HSS approximation, when using pctype ilu (None) 47634f43fa5SPieter Ghysels . -mat_strumpack_reordering <METIS> - Sparsity reducing matrix reordering (choose one of) NATURAL METIS PARMETIS SCOTCH PTSCOTCH RCM (None) 477a2b725a8SWilliam Gropp - -mat_strumpack_iterative_solver <DIRECT> - Select iterative solver from STRUMPACK (choose one of) AUTO DIRECT REFINE PREC_GMRES GMRES PREC_BICGSTAB BICGSTAB (None) 478575704cbSPieter Ghysels 479575704cbSPieter Ghysels Level: beginner 480575704cbSPieter Ghysels 48111a5261eSBarry Smith .seealso: `PCLU`, `PCILU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`, `PCFactorSetMatSolverType()`, `MatSolverType`, `MatGetFactor()` 482575704cbSPieter Ghysels M*/ 4839371c9d4SSatish Balay static PetscErrorCode MatGetFactor_aij_strumpack(Mat A, MatFactorType ftype, Mat *F) { 4847d6ea485SPieter Ghysels Mat B; 4857d6ea485SPieter Ghysels PetscInt M = A->rmap->N, N = A->cmap->N; 48626cc229bSBarry Smith PetscBool verb, flg; 48734f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S; 48834f43fa5SPieter Ghysels STRUMPACK_INTERFACE iface; 4899371c9d4SSatish Balay const STRUMPACK_PRECISION table[2][2][2] = { 4909371c9d4SSatish Balay {{STRUMPACK_FLOATCOMPLEX_64, STRUMPACK_DOUBLECOMPLEX_64}, {STRUMPACK_FLOAT_64, STRUMPACK_DOUBLE_64}}, 4919371c9d4SSatish Balay {{STRUMPACK_FLOATCOMPLEX, STRUMPACK_DOUBLECOMPLEX}, {STRUMPACK_FLOAT, STRUMPACK_DOUBLE} } 4929371c9d4SSatish Balay }; 4934ac6704cSBarry Smith const STRUMPACK_PRECISION prec = table[(sizeof(PetscInt) == 8) ? 0 : 1][(PETSC_SCALAR == PETSC_COMPLEX) ? 0 : 1][(PETSC_REAL == PETSC_FLOAT) ? 0 : 1]; 4947d6ea485SPieter Ghysels 4957d6ea485SPieter Ghysels PetscFunctionBegin; 4967d6ea485SPieter Ghysels /* Create the factorization matrix */ 4979566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B)); 4989566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, M, N)); 4999566063dSJacob Faibussowitsch PetscCall(MatSetType(B, ((PetscObject)A)->type_name)); 5009566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(B, 0, NULL)); 5019566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(B, 0, NULL, 0, NULL)); 50266e17bc3SBarry Smith B->trivialsymbolic = PETSC_TRUE; 503575704cbSPieter Ghysels if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) { 5047d6ea485SPieter Ghysels B->ops->lufactorsymbolic = MatLUFactorSymbolic_STRUMPACK; 505575704cbSPieter Ghysels B->ops->ilufactorsymbolic = MatLUFactorSymbolic_STRUMPACK; 506575704cbSPieter Ghysels } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Factor type not supported"); 5077d6ea485SPieter Ghysels B->ops->view = MatView_STRUMPACK; 5087d6ea485SPieter Ghysels B->ops->destroy = MatDestroy_STRUMPACK; 5097d6ea485SPieter Ghysels B->ops->getdiagonal = MatGetDiagonal_STRUMPACK; 5109566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_aij_strumpack)); 5119566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetReordering_C", MatSTRUMPACKSetReordering_STRUMPACK)); 5129566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetColPerm_C", MatSTRUMPACKSetColPerm_STRUMPACK)); 5139566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetHSSRelTol_C", MatSTRUMPACKSetHSSRelTol_STRUMPACK)); 5149566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetHSSAbsTol_C", MatSTRUMPACKSetHSSAbsTol_STRUMPACK)); 5159566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetHSSMaxRank_C", MatSTRUMPACKSetHSSMaxRank_STRUMPACK)); 5169566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetHSSLeafSize_C", MatSTRUMPACKSetHSSLeafSize_STRUMPACK)); 5179566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetHSSMinSepSize_C", MatSTRUMPACKSetHSSMinSepSize_STRUMPACK)); 518575704cbSPieter Ghysels B->factortype = ftype; 519*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&S)); 5200d08a34dSPieter Ghysels B->spptr = S; 5210d08a34dSPieter Ghysels 5229566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &flg)); 5230d08a34dSPieter Ghysels iface = flg ? STRUMPACK_MT : STRUMPACK_MPI_DIST; 5247d6ea485SPieter Ghysels 52526cc229bSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)B), ((PetscObject)B)->prefix, "STRUMPACK Options", "Mat"); 526575704cbSPieter Ghysels verb = PetscLogPrintInfo ? PETSC_TRUE : PETSC_FALSE; 5279566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_strumpack_verbose", "Print STRUMPACK information", "None", verb, &verb, NULL)); 5287d6ea485SPieter Ghysels 529e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_init", STRUMPACK_init(S, PetscObjectComm((PetscObject)A), prec, iface, 0, NULL, verb)); 53055c022e5SPieter Ghysels 53188382b05SPieter Ghysels if (ftype == MAT_FACTOR_ILU) { 53288382b05SPieter Ghysels /* When enabling HSS compression, the STRUMPACK solver becomes an incomplete */ 53388382b05SPieter Ghysels /* (or approximate) LU factorization. */ 534e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_enable_HSS", STRUMPACK_enable_HSS(*S)); 53588382b05SPieter Ghysels } 536d0609cedSBarry Smith PetscOptionsEnd(); 53755c022e5SPieter Ghysels 53826cc229bSBarry Smith /* set solvertype */ 53926cc229bSBarry Smith PetscCall(PetscFree(B->solvertype)); 54026cc229bSBarry Smith PetscCall(PetscStrallocpy(MATSOLVERSTRUMPACK, &B->solvertype)); 54126cc229bSBarry Smith 5427d6ea485SPieter Ghysels *F = B; 5437d6ea485SPieter Ghysels PetscFunctionReturn(0); 5447d6ea485SPieter Ghysels } 5457d6ea485SPieter Ghysels 5469371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_STRUMPACK(void) { 5477d6ea485SPieter Ghysels PetscFunctionBegin; 5489566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSTRUMPACK, MATMPIAIJ, MAT_FACTOR_LU, MatGetFactor_aij_strumpack)); 5499566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSTRUMPACK, MATSEQAIJ, MAT_FACTOR_LU, MatGetFactor_aij_strumpack)); 5509566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSTRUMPACK, MATMPIAIJ, MAT_FACTOR_ILU, MatGetFactor_aij_strumpack)); 5519566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSTRUMPACK, MATSEQAIJ, MAT_FACTOR_ILU, MatGetFactor_aij_strumpack)); 5527d6ea485SPieter Ghysels PetscFunctionReturn(0); 5537d6ea485SPieter Ghysels } 554