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 5ad0c5e61SPieter Ghysels static PetscErrorCode MatGetDiagonal_STRUMPACK(Mat A,Vec v) 67d6ea485SPieter Ghysels { 77d6ea485SPieter Ghysels PetscFunctionBegin; 87d6ea485SPieter Ghysels SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type: STRUMPACK factor"); 97d6ea485SPieter Ghysels PetscFunctionReturn(0); 107d6ea485SPieter Ghysels } 117d6ea485SPieter Ghysels 12ad0c5e61SPieter Ghysels static PetscErrorCode MatDestroy_STRUMPACK(Mat A) 137d6ea485SPieter Ghysels { 140d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)A->spptr; 157d6ea485SPieter Ghysels PetscBool flg; 167d6ea485SPieter Ghysels 177d6ea485SPieter Ghysels PetscFunctionBegin; 187d6ea485SPieter Ghysels /* Deallocate STRUMPACK storage */ 19*e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_destroy",STRUMPACK_destroy(S)); 209566063dSJacob Faibussowitsch PetscCall(PetscFree(A->spptr)); 219566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg)); 227d6ea485SPieter Ghysels if (flg) { 239566063dSJacob Faibussowitsch PetscCall(MatDestroy_SeqAIJ(A)); 247d6ea485SPieter Ghysels } else { 259566063dSJacob Faibussowitsch PetscCall(MatDestroy_MPIAIJ(A)); 267d6ea485SPieter Ghysels } 27575704cbSPieter Ghysels 28575704cbSPieter Ghysels /* clear composed functions */ 299566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverType_C",NULL)); 309566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetReordering_C",NULL)); 319566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetColPerm_C",NULL)); 329566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSRelTol_C",NULL)); 339566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSAbsTol_C",NULL)); 349566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSMaxRank_C",NULL)); 359566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSLeafSize_C",NULL)); 369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSMinSepSize_C",NULL)); 37575704cbSPieter Ghysels 38575704cbSPieter Ghysels PetscFunctionReturn(0); 39575704cbSPieter Ghysels } 40575704cbSPieter Ghysels 4134f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetReordering_STRUMPACK(Mat F,MatSTRUMPACKReordering reordering) 42575704cbSPieter Ghysels { 430d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 44575704cbSPieter Ghysels 45575704cbSPieter Ghysels PetscFunctionBegin; 46*e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_reordering_method",STRUMPACK_set_reordering_method(*S,(STRUMPACK_REORDERING_STRATEGY)reordering)); 4734f43fa5SPieter Ghysels PetscFunctionReturn(0); 4834f43fa5SPieter Ghysels } 49e5a36eccSStefano Zampini 50542aee0fSPieter Ghysels /*@ 5134f43fa5SPieter Ghysels MatSTRUMPACKSetReordering - Set STRUMPACK fill-reducing reordering 5234f43fa5SPieter Ghysels 5334f43fa5SPieter Ghysels Input Parameters: 5434f43fa5SPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 5534f43fa5SPieter Ghysels - reordering - the code to be used to find the fill-reducing reordering 5634f43fa5SPieter Ghysels Possible values: NATURAL=0 METIS=1 PARMETIS=2 SCOTCH=3 PTSCOTCH=4 RCM=5 5734f43fa5SPieter Ghysels 5834f43fa5SPieter Ghysels Options Database: 5934f43fa5SPieter Ghysels . -mat_strumpack_reordering <METIS> - Sparsity reducing matrix reordering (choose one of) NATURAL METIS PARMETIS SCOTCH PTSCOTCH RCM (None) 6034f43fa5SPieter Ghysels 6134f43fa5SPieter Ghysels Level: beginner 6234f43fa5SPieter Ghysels 6334f43fa5SPieter Ghysels References: 64606c0280SSatish Balay . * - STRUMPACK manual 6534f43fa5SPieter Ghysels 66db781477SPatrick Sanan .seealso: `MatGetFactor()` 67542aee0fSPieter Ghysels @*/ 6834f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetReordering(Mat F,MatSTRUMPACKReordering reordering) 6934f43fa5SPieter Ghysels { 7034f43fa5SPieter Ghysels PetscFunctionBegin; 7134f43fa5SPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 7234f43fa5SPieter Ghysels PetscValidLogicalCollectiveEnum(F,reordering,2); 73cac4c232SBarry Smith PetscTryMethod(F,"MatSTRUMPACKSetReordering_C",(Mat,MatSTRUMPACKReordering),(F,reordering)); 74575704cbSPieter Ghysels PetscFunctionReturn(0); 75575704cbSPieter Ghysels } 76575704cbSPieter Ghysels 7734f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetColPerm_STRUMPACK(Mat F,PetscBool cperm) 7834f43fa5SPieter Ghysels { 7934f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 8034f43fa5SPieter Ghysels 8134f43fa5SPieter Ghysels PetscFunctionBegin; 82*e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_set_mc64job",STRUMPACK_set_mc64job(*S,cperm ? 5 : 0)); 8334f43fa5SPieter Ghysels PetscFunctionReturn(0); 8434f43fa5SPieter Ghysels } 85e5a36eccSStefano Zampini 86575704cbSPieter Ghysels /*@ 8734f43fa5SPieter Ghysels MatSTRUMPACKSetColPerm - Set whether STRUMPACK should try to permute the columns of the matrix in order to get a nonzero diagonal 88147403d9SBarry Smith 89575704cbSPieter Ghysels Logically Collective on Mat 90575704cbSPieter Ghysels 91575704cbSPieter Ghysels Input Parameters: 92575704cbSPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 93147403d9SBarry Smith - cperm - PETSC_TRUE to permute (internally) the columns of the matrix 94575704cbSPieter Ghysels 95575704cbSPieter Ghysels Options Database: 96147403d9SBarry Smith . -mat_strumpack_colperm <cperm> - true to use the permutation 97575704cbSPieter Ghysels 98575704cbSPieter Ghysels Level: beginner 99575704cbSPieter Ghysels 100575704cbSPieter Ghysels References: 101606c0280SSatish Balay . * - STRUMPACK manual 102575704cbSPieter Ghysels 103db781477SPatrick Sanan .seealso: `MatGetFactor()` 104575704cbSPieter Ghysels @*/ 10534f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetColPerm(Mat F,PetscBool cperm) 106575704cbSPieter Ghysels { 107575704cbSPieter Ghysels PetscFunctionBegin; 108575704cbSPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 10934f43fa5SPieter Ghysels PetscValidLogicalCollectiveBool(F,cperm,2); 110cac4c232SBarry Smith PetscTryMethod(F,"MatSTRUMPACKSetColPerm_C",(Mat,PetscBool),(F,cperm)); 111575704cbSPieter Ghysels PetscFunctionReturn(0); 112575704cbSPieter Ghysels } 113575704cbSPieter Ghysels 11434f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSRelTol_STRUMPACK(Mat F,PetscReal rtol) 11534f43fa5SPieter Ghysels { 11634f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 11734f43fa5SPieter Ghysels 11834f43fa5SPieter Ghysels PetscFunctionBegin; 119*e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_set_HSS_rel_tol", STRUMPACK_set_HSS_rel_tol(*S,rtol)); 12034f43fa5SPieter Ghysels PetscFunctionReturn(0); 12134f43fa5SPieter Ghysels } 122e5a36eccSStefano Zampini 12334f43fa5SPieter Ghysels /*@ 12434f43fa5SPieter Ghysels MatSTRUMPACKSetHSSRelTol - Set STRUMPACK relative tolerance for HSS compression 125147403d9SBarry Smith 12634f43fa5SPieter Ghysels Logically Collective on Mat 12734f43fa5SPieter Ghysels 12834f43fa5SPieter Ghysels Input Parameters: 12934f43fa5SPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 13034f43fa5SPieter Ghysels - rtol - relative compression tolerance 13134f43fa5SPieter Ghysels 13234f43fa5SPieter Ghysels Options Database: 13334f43fa5SPieter Ghysels . -mat_strumpack_hss_rel_tol <1e-2> - Relative compression tolerance (None) 13434f43fa5SPieter Ghysels 13534f43fa5SPieter Ghysels Level: beginner 13634f43fa5SPieter Ghysels 13734f43fa5SPieter Ghysels References: 138606c0280SSatish Balay . * - STRUMPACK manual 13934f43fa5SPieter Ghysels 140db781477SPatrick Sanan .seealso: `MatGetFactor()` 14134f43fa5SPieter Ghysels @*/ 14234f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSRelTol(Mat F,PetscReal rtol) 14334f43fa5SPieter Ghysels { 14434f43fa5SPieter Ghysels PetscFunctionBegin; 14534f43fa5SPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 14634f43fa5SPieter Ghysels PetscValidLogicalCollectiveReal(F,rtol,2); 147cac4c232SBarry Smith PetscTryMethod(F,"MatSTRUMPACKSetHSSRelTol_C",(Mat,PetscReal),(F,rtol)); 14834f43fa5SPieter Ghysels PetscFunctionReturn(0); 14934f43fa5SPieter Ghysels } 15034f43fa5SPieter Ghysels 15134f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSAbsTol_STRUMPACK(Mat F,PetscReal atol) 15234f43fa5SPieter Ghysels { 15334f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 15434f43fa5SPieter Ghysels 15534f43fa5SPieter Ghysels PetscFunctionBegin; 156*e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_set_HSS_abs_tol", STRUMPACK_set_HSS_abs_tol(*S,atol)); 15734f43fa5SPieter Ghysels PetscFunctionReturn(0); 15834f43fa5SPieter Ghysels } 159e5a36eccSStefano Zampini 16034f43fa5SPieter Ghysels /*@ 16134f43fa5SPieter Ghysels MatSTRUMPACKSetHSSAbsTol - Set STRUMPACK absolute tolerance for HSS compression 162147403d9SBarry Smith 16334f43fa5SPieter Ghysels Logically Collective on Mat 16434f43fa5SPieter Ghysels 16534f43fa5SPieter Ghysels Input Parameters: 16634f43fa5SPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 16734f43fa5SPieter Ghysels - atol - absolute compression tolerance 16834f43fa5SPieter Ghysels 16934f43fa5SPieter Ghysels Options Database: 17034f43fa5SPieter Ghysels . -mat_strumpack_hss_abs_tol <1e-8> - Absolute compression tolerance (None) 17134f43fa5SPieter Ghysels 17234f43fa5SPieter Ghysels Level: beginner 17334f43fa5SPieter Ghysels 17434f43fa5SPieter Ghysels References: 175606c0280SSatish Balay . * - STRUMPACK manual 17634f43fa5SPieter Ghysels 177db781477SPatrick Sanan .seealso: `MatGetFactor()` 17834f43fa5SPieter Ghysels @*/ 17934f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSAbsTol(Mat F,PetscReal atol) 18034f43fa5SPieter Ghysels { 18134f43fa5SPieter Ghysels PetscFunctionBegin; 18234f43fa5SPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 18334f43fa5SPieter Ghysels PetscValidLogicalCollectiveReal(F,atol,2); 184cac4c232SBarry Smith PetscTryMethod(F,"MatSTRUMPACKSetHSSAbsTol_C",(Mat,PetscReal),(F,atol)); 18534f43fa5SPieter Ghysels PetscFunctionReturn(0); 18634f43fa5SPieter Ghysels } 18734f43fa5SPieter Ghysels 18834f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSMaxRank_STRUMPACK(Mat F,PetscInt hssmaxrank) 18934f43fa5SPieter Ghysels { 19034f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 19134f43fa5SPieter Ghysels 19234f43fa5SPieter Ghysels PetscFunctionBegin; 193*e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_set_HSS_max_rank", STRUMPACK_set_HSS_max_rank(*S,hssmaxrank)); 19434f43fa5SPieter Ghysels PetscFunctionReturn(0); 19534f43fa5SPieter Ghysels } 196e5a36eccSStefano Zampini 19734f43fa5SPieter Ghysels /*@ 19834f43fa5SPieter Ghysels MatSTRUMPACKSetHSSMaxRank - Set STRUMPACK maximum HSS rank 199147403d9SBarry Smith 20034f43fa5SPieter Ghysels Logically Collective on Mat 20134f43fa5SPieter Ghysels 20234f43fa5SPieter Ghysels Input Parameters: 20334f43fa5SPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 20434f43fa5SPieter Ghysels - hssmaxrank - maximum rank used in low-rank approximation 20534f43fa5SPieter Ghysels 20634f43fa5SPieter Ghysels Options Database: 20734f43fa5SPieter Ghysels . -mat_strumpack_max_rank - Maximum rank in HSS compression, when using pctype ilu (None) 20834f43fa5SPieter Ghysels 20934f43fa5SPieter Ghysels Level: beginner 21034f43fa5SPieter Ghysels 21134f43fa5SPieter Ghysels References: 212606c0280SSatish Balay . * - STRUMPACK manual 21334f43fa5SPieter Ghysels 214db781477SPatrick Sanan .seealso: `MatGetFactor()` 21534f43fa5SPieter Ghysels @*/ 21634f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSMaxRank(Mat F,PetscInt hssmaxrank) 21734f43fa5SPieter Ghysels { 21834f43fa5SPieter Ghysels PetscFunctionBegin; 21934f43fa5SPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 22034f43fa5SPieter Ghysels PetscValidLogicalCollectiveInt(F,hssmaxrank,2); 221cac4c232SBarry Smith PetscTryMethod(F,"MatSTRUMPACKSetHSSMaxRank_C",(Mat,PetscInt),(F,hssmaxrank)); 22234f43fa5SPieter Ghysels PetscFunctionReturn(0); 22334f43fa5SPieter Ghysels } 22434f43fa5SPieter Ghysels 225a36bf211SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSLeafSize_STRUMPACK(Mat F,PetscInt leaf_size) 226a36bf211SPieter Ghysels { 227a36bf211SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 228a36bf211SPieter Ghysels 229a36bf211SPieter Ghysels PetscFunctionBegin; 230*e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_set_HSS_leaf_size", STRUMPACK_set_HSS_leaf_size(*S,leaf_size)); 231a36bf211SPieter Ghysels PetscFunctionReturn(0); 232a36bf211SPieter Ghysels } 233e5a36eccSStefano Zampini 234a36bf211SPieter Ghysels /*@ 235a36bf211SPieter Ghysels MatSTRUMPACKSetHSSLeafSize - Set STRUMPACK HSS leaf size 236147403d9SBarry Smith 237a36bf211SPieter Ghysels Logically Collective on Mat 238a36bf211SPieter Ghysels 239a36bf211SPieter Ghysels Input Parameters: 240a36bf211SPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 241a36bf211SPieter Ghysels - leaf_size - Size of diagonal blocks in HSS approximation 242a36bf211SPieter Ghysels 243a36bf211SPieter Ghysels Options Database: 244a36bf211SPieter Ghysels . -mat_strumpack_leaf_size - Size of diagonal blocks in HSS approximation, when using pctype ilu (None) 245a36bf211SPieter Ghysels 246a36bf211SPieter Ghysels Level: beginner 247a36bf211SPieter Ghysels 248a36bf211SPieter Ghysels References: 249606c0280SSatish Balay . * - STRUMPACK manual 250a36bf211SPieter Ghysels 251db781477SPatrick Sanan .seealso: `MatGetFactor()` 252a36bf211SPieter Ghysels @*/ 253a36bf211SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSLeafSize(Mat F,PetscInt leaf_size) 254a36bf211SPieter Ghysels { 255a36bf211SPieter Ghysels PetscFunctionBegin; 256a36bf211SPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 257a36bf211SPieter Ghysels PetscValidLogicalCollectiveInt(F,leaf_size,2); 258cac4c232SBarry Smith PetscTryMethod(F,"MatSTRUMPACKSetHSSLeafSize_C",(Mat,PetscInt),(F,leaf_size)); 259a36bf211SPieter Ghysels PetscFunctionReturn(0); 260a36bf211SPieter Ghysels } 261a36bf211SPieter Ghysels 26234f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSMinSepSize_STRUMPACK(Mat F,PetscInt hssminsize) 26334f43fa5SPieter Ghysels { 26434f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 26534f43fa5SPieter Ghysels 26634f43fa5SPieter Ghysels PetscFunctionBegin; 267*e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_set_HSS_min_sep_size", STRUMPACK_set_HSS_min_sep_size(*S,hssminsize)); 26834f43fa5SPieter Ghysels PetscFunctionReturn(0); 26934f43fa5SPieter Ghysels } 270e5a36eccSStefano Zampini 271291d0ed5SPieter Ghysels /*@ 272291d0ed5SPieter Ghysels MatSTRUMPACKSetHSSMinSepSize - Set STRUMPACK minimum separator size for low-rank approximation 273147403d9SBarry Smith 274291d0ed5SPieter Ghysels Logically Collective on Mat 275291d0ed5SPieter Ghysels 276291d0ed5SPieter Ghysels Input Parameters: 277291d0ed5SPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 278291d0ed5SPieter Ghysels - hssminsize - minimum dense matrix size for low-rank approximation 279291d0ed5SPieter Ghysels 280291d0ed5SPieter Ghysels Options Database: 281147403d9SBarry Smith . -mat_strumpack_hss_min_sep_size <hssminsize> - set the minimum separator size 282291d0ed5SPieter Ghysels 283291d0ed5SPieter Ghysels Level: beginner 284291d0ed5SPieter Ghysels 285291d0ed5SPieter Ghysels References: 286606c0280SSatish Balay . * - STRUMPACK manual 287291d0ed5SPieter Ghysels 288db781477SPatrick Sanan .seealso: `MatGetFactor()` 289291d0ed5SPieter Ghysels @*/ 290291d0ed5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat F,PetscInt hssminsize) 291291d0ed5SPieter Ghysels { 292291d0ed5SPieter Ghysels PetscFunctionBegin; 293291d0ed5SPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 294291d0ed5SPieter Ghysels PetscValidLogicalCollectiveInt(F,hssminsize,2); 295cac4c232SBarry Smith PetscTryMethod(F,"MatSTRUMPACKSetHSSMinSepSize_C",(Mat,PetscInt),(F,hssminsize)); 296291d0ed5SPieter Ghysels PetscFunctionReturn(0); 297291d0ed5SPieter Ghysels } 298291d0ed5SPieter Ghysels 299ad0c5e61SPieter Ghysels static PetscErrorCode MatSolve_STRUMPACK(Mat A,Vec b_mpi,Vec x) 3007d6ea485SPieter Ghysels { 3010d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)A->spptr; 3027d6ea485SPieter Ghysels STRUMPACK_RETURN_CODE sp_err; 3037d6ea485SPieter Ghysels const PetscScalar *bptr; 3047d6ea485SPieter Ghysels PetscScalar *xptr; 3057d6ea485SPieter Ghysels 3067d6ea485SPieter Ghysels PetscFunctionBegin; 3079566063dSJacob Faibussowitsch PetscCall(VecGetArray(x,&xptr)); 3089566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(b_mpi,&bptr)); 3090d08a34dSPieter Ghysels 310*e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_solve",sp_err = STRUMPACK_solve(*S,(PetscScalar*)bptr,xptr,0)); 3110d08a34dSPieter Ghysels switch (sp_err) { 3120d08a34dSPieter Ghysels case STRUMPACK_SUCCESS: break; 3130d08a34dSPieter Ghysels case STRUMPACK_MATRIX_NOT_SET: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); break; } 3140d08a34dSPieter Ghysels case STRUMPACK_REORDERING_ERROR: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); break; } 3150d08a34dSPieter Ghysels default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: solve failed"); 3167d6ea485SPieter Ghysels } 3179566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x,&xptr)); 3189566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(b_mpi,&bptr)); 3197d6ea485SPieter Ghysels PetscFunctionReturn(0); 3207d6ea485SPieter Ghysels } 3217d6ea485SPieter Ghysels 322ad0c5e61SPieter Ghysels static PetscErrorCode MatMatSolve_STRUMPACK(Mat A,Mat B_mpi,Mat X) 3237d6ea485SPieter Ghysels { 3247d6ea485SPieter Ghysels PetscBool flg; 3257d6ea485SPieter Ghysels 3267d6ea485SPieter Ghysels PetscFunctionBegin; 3279566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)B_mpi,&flg,MATSEQDENSE,MATMPIDENSE,NULL)); 3285f80ce2aSJacob Faibussowitsch PetscCheck(flg,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 3299566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL)); 3305f80ce2aSJacob Faibussowitsch PetscCheck(flg,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); 3317d6ea485SPieter Ghysels SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_STRUMPACK() is not implemented yet"); 3327d6ea485SPieter Ghysels PetscFunctionReturn(0); 3337d6ea485SPieter Ghysels } 3347d6ea485SPieter Ghysels 335860c79edSBarry Smith static PetscErrorCode MatView_Info_STRUMPACK(Mat A,PetscViewer viewer) 336ad0c5e61SPieter Ghysels { 337ad0c5e61SPieter Ghysels PetscFunctionBegin; 338ad0c5e61SPieter Ghysels /* check if matrix is strumpack type */ 339ad0c5e61SPieter Ghysels if (A->ops->solve != MatSolve_STRUMPACK) PetscFunctionReturn(0); 3409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"STRUMPACK sparse solver!\n")); 341ad0c5e61SPieter Ghysels PetscFunctionReturn(0); 342ad0c5e61SPieter Ghysels } 343ad0c5e61SPieter Ghysels 344ad0c5e61SPieter Ghysels static PetscErrorCode MatView_STRUMPACK(Mat A,PetscViewer viewer) 345ad0c5e61SPieter Ghysels { 346ad0c5e61SPieter Ghysels PetscBool iascii; 347ad0c5e61SPieter Ghysels PetscViewerFormat format; 348ad0c5e61SPieter Ghysels 349ad0c5e61SPieter Ghysels PetscFunctionBegin; 3509566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 351ad0c5e61SPieter Ghysels if (iascii) { 3529566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer,&format)); 3539566063dSJacob Faibussowitsch if (format == PETSC_VIEWER_ASCII_INFO) PetscCall(MatView_Info_STRUMPACK(A,viewer)); 354ad0c5e61SPieter Ghysels } 355ad0c5e61SPieter Ghysels PetscFunctionReturn(0); 356ad0c5e61SPieter Ghysels } 3577d6ea485SPieter Ghysels 358ad0c5e61SPieter Ghysels static PetscErrorCode MatLUFactorNumeric_STRUMPACK(Mat F,Mat A,const MatFactorInfo *info) 3597d6ea485SPieter Ghysels { 3600d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 3617d6ea485SPieter Ghysels STRUMPACK_RETURN_CODE sp_err; 3627d6ea485SPieter Ghysels Mat_SeqAIJ *A_d,*A_o; 3637d6ea485SPieter Ghysels Mat_MPIAIJ *mat; 3640d08a34dSPieter Ghysels PetscInt M=A->rmap->N,m=A->rmap->n; 3657d6ea485SPieter Ghysels PetscBool flg; 3667d6ea485SPieter Ghysels 3677d6ea485SPieter Ghysels PetscFunctionBegin; 3689566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg)); 3697d6ea485SPieter Ghysels if (flg) { /* A is MATMPIAIJ */ 3707d6ea485SPieter Ghysels mat = (Mat_MPIAIJ*)A->data; 3717d6ea485SPieter Ghysels A_d = (Mat_SeqAIJ*)(mat->A)->data; 3727d6ea485SPieter Ghysels A_o = (Mat_SeqAIJ*)(mat->B)->data; 373*e77caa6dSBarry 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)); 3747d6ea485SPieter Ghysels } else { /* A is MATSEQAIJ */ 3757d6ea485SPieter Ghysels A_d = (Mat_SeqAIJ*)A->data; 376*e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_set_csr_matrix",STRUMPACK_set_csr_matrix(*S,&M,A_d->i,A_d->j,A_d->a,0)); 3777d6ea485SPieter Ghysels } 3787d6ea485SPieter Ghysels 3797d6ea485SPieter Ghysels /* Reorder and Factor the matrix. */ 3807d6ea485SPieter Ghysels /* TODO figure out how to avoid reorder if the matrix values changed, but the pattern remains the same. */ 381*e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_reorder",sp_err = STRUMPACK_reorder(*S)); 382*e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_factor",sp_err = STRUMPACK_factor(*S)); 3830d08a34dSPieter Ghysels switch (sp_err) { 3840d08a34dSPieter Ghysels case STRUMPACK_SUCCESS: break; 3850d08a34dSPieter Ghysels case STRUMPACK_MATRIX_NOT_SET: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); break; } 3860d08a34dSPieter Ghysels case STRUMPACK_REORDERING_ERROR: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); break; } 3870d08a34dSPieter Ghysels default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: factorization failed"); 3887d6ea485SPieter Ghysels } 389cb250fa3SPieter Ghysels F->assembled = PETSC_TRUE; 390cb250fa3SPieter Ghysels F->preallocated = PETSC_TRUE; 3917d6ea485SPieter Ghysels PetscFunctionReturn(0); 3927d6ea485SPieter Ghysels } 3937d6ea485SPieter Ghysels 394ad0c5e61SPieter Ghysels static PetscErrorCode MatLUFactorSymbolic_STRUMPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 3957d6ea485SPieter Ghysels { 39626cc229bSBarry Smith STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 39726cc229bSBarry Smith PetscBool flg,set; 39826cc229bSBarry Smith PetscReal ctol; 39926cc229bSBarry Smith PetscInt hssminsize,max_rank,leaf_size; 40026cc229bSBarry Smith STRUMPACK_REORDERING_STRATEGY ndcurrent,ndvalue; 40126cc229bSBarry Smith STRUMPACK_KRYLOV_SOLVER itcurrent,itsolver; 40226cc229bSBarry Smith const char *const STRUMPACKNDTypes[] = {"NATURAL","METIS","PARMETIS","SCOTCH","PTSCOTCH","RCM","STRUMPACKNDTypes","",0}; 40326cc229bSBarry Smith const char *const SolverTypes[] = {"AUTO","NONE","REFINE","PREC_GMRES","GMRES","PREC_BICGSTAB","BICGSTAB","SolverTypes","",0}; 40426cc229bSBarry Smith 4057d6ea485SPieter Ghysels PetscFunctionBegin; 40626cc229bSBarry Smith /* Set options to F */ 40726cc229bSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)F),((PetscObject)F)->prefix,"STRUMPACK Options","Mat"); 408*e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_HSS_rel_tol",ctol = (PetscReal)STRUMPACK_HSS_rel_tol(*S)); 40926cc229bSBarry Smith PetscCall(PetscOptionsReal("-mat_strumpack_hss_rel_tol","Relative compression tolerance","None",ctol,&ctol,&set)); 410*e77caa6dSBarry Smith if (set) PetscStackCallExternalVoid("STRUMPACK_set_HSS_rel_tol",STRUMPACK_set_HSS_rel_tol(*S,(double)ctol)); 41126cc229bSBarry Smith 412*e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_HSS_abs_tol",ctol = (PetscReal)STRUMPACK_HSS_abs_tol(*S)); 41326cc229bSBarry Smith PetscCall(PetscOptionsReal("-mat_strumpack_hss_abs_tol","Absolute compression tolerance","None",ctol,&ctol,&set)); 414*e77caa6dSBarry Smith if (set) PetscStackCallExternalVoid("STRUMPACK_set_HSS_abs_tol",STRUMPACK_set_HSS_abs_tol(*S,(double)ctol)); 41526cc229bSBarry Smith 416*e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_mc64job",flg = (STRUMPACK_mc64job(*S) == 0) ? PETSC_FALSE : PETSC_TRUE); 41726cc229bSBarry Smith PetscCall(PetscOptionsBool("-mat_strumpack_colperm","Find a col perm to get nonzero diagonal","None",flg,&flg,&set)); 418*e77caa6dSBarry Smith if (set) PetscStackCallExternalVoid("STRUMPACK_set_mc64job",STRUMPACK_set_mc64job(*S,flg ? 5 : 0)); 41926cc229bSBarry Smith 420*e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_HSS_min_sep_size",hssminsize = (PetscInt)STRUMPACK_HSS_min_sep_size(*S)); 42126cc229bSBarry Smith PetscCall(PetscOptionsInt("-mat_strumpack_hss_min_sep_size","Minimum size of separator for HSS compression","None",hssminsize,&hssminsize,&set)); 422*e77caa6dSBarry Smith if (set) PetscStackCallExternalVoid("STRUMPACK_set_HSS_min_sep_size",STRUMPACK_set_HSS_min_sep_size(*S,(int)hssminsize)); 42326cc229bSBarry Smith 424*e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_HSS_max_rank",max_rank = (PetscInt)STRUMPACK_HSS_max_rank(*S)); 42526cc229bSBarry Smith PetscCall(PetscOptionsInt("-mat_strumpack_max_rank","Maximum rank in HSS approximation","None",max_rank,&max_rank,&set)); 426*e77caa6dSBarry Smith if (set) PetscStackCallExternalVoid("STRUMPACK_set_HSS_max_rank",STRUMPACK_set_HSS_max_rank(*S,(int)max_rank)); 42726cc229bSBarry Smith 428*e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_HSS_leaf_size",leaf_size = (PetscInt)STRUMPACK_HSS_leaf_size(*S)); 42926cc229bSBarry Smith PetscCall(PetscOptionsInt("-mat_strumpack_leaf_size","Size of diagonal blocks in HSS approximation","None",leaf_size,&leaf_size,&set)); 430*e77caa6dSBarry Smith if (set) PetscStackCallExternalVoid("STRUMPACK_set_HSS_leaf_size",STRUMPACK_set_HSS_leaf_size(*S,(int)leaf_size)); 43126cc229bSBarry Smith 432*e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_reordering_method",ndcurrent = STRUMPACK_reordering_method(*S)); 43326cc229bSBarry Smith PetscCall(PetscOptionsEnum("-mat_strumpack_reordering","Sparsity reducing matrix reordering","None",STRUMPACKNDTypes,(PetscEnum)ndcurrent,(PetscEnum*)&ndvalue,&set)); 434*e77caa6dSBarry Smith if (set) PetscStackCallExternalVoid("STRUMPACK_set_reordering_method",STRUMPACK_set_reordering_method(*S,ndvalue)); 43526cc229bSBarry Smith 43626cc229bSBarry Smith /* Disable the outer iterative solver from STRUMPACK. */ 43726cc229bSBarry Smith /* When STRUMPACK is used as a direct solver, it will by default do iterative refinement. */ 43826cc229bSBarry Smith /* When STRUMPACK is used as an approximate factorization preconditioner (by enabling */ 43926cc229bSBarry Smith /* low-rank compression), it will use it's own preconditioned GMRES. Here we can disable */ 44026cc229bSBarry Smith /* the outer iterative solver, as PETSc uses STRUMPACK from within a KSP. */ 441*e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_set_Krylov_solver", STRUMPACK_set_Krylov_solver(*S, STRUMPACK_DIRECT)); 44226cc229bSBarry Smith 443*e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_Krylov_solver",itcurrent = STRUMPACK_Krylov_solver(*S)); 44426cc229bSBarry Smith PetscCall(PetscOptionsEnum("-mat_strumpack_iterative_solver","Select iterative solver from STRUMPACK","None",SolverTypes,(PetscEnum)itcurrent,(PetscEnum*)&itsolver,&set)); 445*e77caa6dSBarry Smith if (set) PetscStackCallExternalVoid("STRUMPACK_set_Krylov_solver",STRUMPACK_set_Krylov_solver(*S,itsolver)); 44626cc229bSBarry Smith PetscOptionsEnd(); 44726cc229bSBarry Smith 4487d6ea485SPieter Ghysels F->ops->lufactornumeric = MatLUFactorNumeric_STRUMPACK; 4497d6ea485SPieter Ghysels F->ops->solve = MatSolve_STRUMPACK; 4507d6ea485SPieter Ghysels F->ops->matsolve = MatMatSolve_STRUMPACK; 4517d6ea485SPieter Ghysels PetscFunctionReturn(0); 4527d6ea485SPieter Ghysels } 4537d6ea485SPieter Ghysels 454ea799195SBarry Smith static PetscErrorCode MatFactorGetSolverType_aij_strumpack(Mat A,MatSolverType *type) 4557d6ea485SPieter Ghysels { 4567d6ea485SPieter Ghysels PetscFunctionBegin; 4577d6ea485SPieter Ghysels *type = MATSOLVERSTRUMPACK; 4587d6ea485SPieter Ghysels PetscFunctionReturn(0); 4597d6ea485SPieter Ghysels } 4607d6ea485SPieter Ghysels 461575704cbSPieter Ghysels /*MC 462575704cbSPieter Ghysels MATSOLVERSSTRUMPACK = "strumpack" - A solver package providing a direct sparse solver (PCLU) 463575704cbSPieter Ghysels and a preconditioner (PCILU) using low-rank compression via the external package STRUMPACK. 464575704cbSPieter Ghysels 465575704cbSPieter Ghysels Consult the STRUMPACK-sparse manual for more info. 466575704cbSPieter Ghysels 467575704cbSPieter Ghysels Use 468575704cbSPieter Ghysels ./configure --download-strumpack 469575704cbSPieter Ghysels to have PETSc installed with STRUMPACK 470575704cbSPieter Ghysels 471575704cbSPieter Ghysels Use 4723ca39a21SBarry Smith -pc_type lu -pc_factor_mat_solver_type strumpack 473575704cbSPieter Ghysels to use this as an exact (direct) solver, use 4743ca39a21SBarry Smith -pc_type ilu -pc_factor_mat_solver_type strumpack 475575704cbSPieter Ghysels to enable low-rank compression (i.e, use as a preconditioner). 476575704cbSPieter Ghysels 477575704cbSPieter Ghysels Works with AIJ matrices 478575704cbSPieter Ghysels 479575704cbSPieter Ghysels Options Database Keys: 48067b8a455SSatish Balay + -mat_strumpack_verbose - verbose info 48134f43fa5SPieter Ghysels . -mat_strumpack_hss_rel_tol <1e-2> - Relative compression tolerance (None) 48234f43fa5SPieter Ghysels . -mat_strumpack_hss_abs_tol <1e-8> - Absolute compression tolerance (None) 483575704cbSPieter Ghysels . -mat_strumpack_colperm <TRUE> - Permute matrix to make diagonal nonzeros (None) 48434f43fa5SPieter Ghysels . -mat_strumpack_hss_min_sep_size <256> - Minimum size of separator for HSS compression (None) 48534f43fa5SPieter Ghysels . -mat_strumpack_max_rank - Maximum rank in HSS compression, when using pctype ilu (None) 486a36bf211SPieter Ghysels . -mat_strumpack_leaf_size - Size of diagonal blocks in HSS approximation, when using pctype ilu (None) 48734f43fa5SPieter Ghysels . -mat_strumpack_reordering <METIS> - Sparsity reducing matrix reordering (choose one of) NATURAL METIS PARMETIS SCOTCH PTSCOTCH RCM (None) 488a2b725a8SWilliam Gropp - -mat_strumpack_iterative_solver <DIRECT> - Select iterative solver from STRUMPACK (choose one of) AUTO DIRECT REFINE PREC_GMRES GMRES PREC_BICGSTAB BICGSTAB (None) 489575704cbSPieter Ghysels 490575704cbSPieter Ghysels Level: beginner 491575704cbSPieter Ghysels 492db781477SPatrick Sanan .seealso: `PCLU`, `PCILU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`, `PCFactorSetMatSolverType()`, `MatSolverType` 493575704cbSPieter Ghysels M*/ 494ad0c5e61SPieter Ghysels static PetscErrorCode MatGetFactor_aij_strumpack(Mat A,MatFactorType ftype,Mat *F) 4957d6ea485SPieter Ghysels { 4967d6ea485SPieter Ghysels Mat B; 4977d6ea485SPieter Ghysels PetscInt M=A->rmap->N,N=A->cmap->N; 49826cc229bSBarry Smith PetscBool verb,flg; 49934f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S; 50034f43fa5SPieter Ghysels STRUMPACK_INTERFACE iface; 50134f43fa5SPieter Ghysels const STRUMPACK_PRECISION table[2][2][2] = 50234f43fa5SPieter Ghysels {{{STRUMPACK_FLOATCOMPLEX_64, STRUMPACK_DOUBLECOMPLEX_64}, 5037d6ea485SPieter Ghysels {STRUMPACK_FLOAT_64, STRUMPACK_DOUBLE_64}}, 5047d6ea485SPieter Ghysels {{STRUMPACK_FLOATCOMPLEX, STRUMPACK_DOUBLECOMPLEX}, 5057d6ea485SPieter Ghysels {STRUMPACK_FLOAT, STRUMPACK_DOUBLE}}}; 5064ac6704cSBarry Smith const STRUMPACK_PRECISION prec = table[(sizeof(PetscInt)==8)?0:1][(PETSC_SCALAR==PETSC_COMPLEX)?0:1][(PETSC_REAL==PETSC_FLOAT)?0:1]; 5077d6ea485SPieter Ghysels 5087d6ea485SPieter Ghysels PetscFunctionBegin; 5097d6ea485SPieter Ghysels /* Create the factorization matrix */ 5109566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&B)); 5119566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B,A->rmap->n,A->cmap->n,M,N)); 5129566063dSJacob Faibussowitsch PetscCall(MatSetType(B,((PetscObject)A)->type_name)); 5139566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(B,0,NULL)); 5149566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(B,0,NULL,0,NULL)); 51566e17bc3SBarry Smith B->trivialsymbolic = PETSC_TRUE; 516575704cbSPieter Ghysels if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) { 5177d6ea485SPieter Ghysels B->ops->lufactorsymbolic = MatLUFactorSymbolic_STRUMPACK; 518575704cbSPieter Ghysels B->ops->ilufactorsymbolic = MatLUFactorSymbolic_STRUMPACK; 519575704cbSPieter Ghysels } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 5207d6ea485SPieter Ghysels B->ops->view = MatView_STRUMPACK; 5217d6ea485SPieter Ghysels B->ops->destroy = MatDestroy_STRUMPACK; 5227d6ea485SPieter Ghysels B->ops->getdiagonal = MatGetDiagonal_STRUMPACK; 5239566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_aij_strumpack)); 5249566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetReordering_C",MatSTRUMPACKSetReordering_STRUMPACK)); 5259566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetColPerm_C",MatSTRUMPACKSetColPerm_STRUMPACK)); 5269566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSRelTol_C",MatSTRUMPACKSetHSSRelTol_STRUMPACK)); 5279566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSAbsTol_C",MatSTRUMPACKSetHSSAbsTol_STRUMPACK)); 5289566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSMaxRank_C",MatSTRUMPACKSetHSSMaxRank_STRUMPACK)); 5299566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSLeafSize_C",MatSTRUMPACKSetHSSLeafSize_STRUMPACK)); 5309566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSMinSepSize_C",MatSTRUMPACKSetHSSMinSepSize_STRUMPACK)); 531575704cbSPieter Ghysels B->factortype = ftype; 5329566063dSJacob Faibussowitsch PetscCall(PetscNewLog(B,&S)); 5330d08a34dSPieter Ghysels B->spptr = S; 5340d08a34dSPieter Ghysels 5359566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg)); 5360d08a34dSPieter Ghysels iface = flg ? STRUMPACK_MT : STRUMPACK_MPI_DIST; 5377d6ea485SPieter Ghysels 53826cc229bSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)B),((PetscObject)B)->prefix,"STRUMPACK Options","Mat"); 539575704cbSPieter Ghysels verb = PetscLogPrintInfo ? PETSC_TRUE : PETSC_FALSE; 5409566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_strumpack_verbose","Print STRUMPACK information","None",verb,&verb,NULL)); 5417d6ea485SPieter Ghysels 542*e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_init",STRUMPACK_init(S,PetscObjectComm((PetscObject)A),prec,iface,0,NULL,verb)); 54355c022e5SPieter Ghysels 54488382b05SPieter Ghysels if (ftype == MAT_FACTOR_ILU) { 54588382b05SPieter Ghysels /* When enabling HSS compression, the STRUMPACK solver becomes an incomplete */ 54688382b05SPieter Ghysels /* (or approximate) LU factorization. */ 547*e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_enable_HSS",STRUMPACK_enable_HSS(*S)); 54888382b05SPieter Ghysels } 549d0609cedSBarry Smith PetscOptionsEnd(); 55055c022e5SPieter Ghysels 55126cc229bSBarry Smith /* set solvertype */ 55226cc229bSBarry Smith PetscCall(PetscFree(B->solvertype)); 55326cc229bSBarry Smith PetscCall(PetscStrallocpy(MATSOLVERSTRUMPACK,&B->solvertype)); 55426cc229bSBarry Smith 5557d6ea485SPieter Ghysels *F = B; 5567d6ea485SPieter Ghysels PetscFunctionReturn(0); 5577d6ea485SPieter Ghysels } 5587d6ea485SPieter Ghysels 5593ca39a21SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_STRUMPACK(void) 5607d6ea485SPieter Ghysels { 5617d6ea485SPieter Ghysels PetscFunctionBegin; 5629566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack)); 5639566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack)); 5649566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_ILU,MatGetFactor_aij_strumpack)); 5659566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_ILU,MatGetFactor_aij_strumpack)); 5667d6ea485SPieter Ghysels PetscFunctionReturn(0); 5677d6ea485SPieter Ghysels } 568