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 5291d0ed5SPieter Ghysels #undef __FUNCT__ 6291d0ed5SPieter Ghysels #define __FUNCT__ "MatGetDiagonal_STRUMPACK" 7ad0c5e61SPieter Ghysels static PetscErrorCode MatGetDiagonal_STRUMPACK(Mat A,Vec v) 87d6ea485SPieter Ghysels { 97d6ea485SPieter Ghysels PetscFunctionBegin; 107d6ea485SPieter Ghysels SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type: STRUMPACK factor"); 117d6ea485SPieter Ghysels PetscFunctionReturn(0); 127d6ea485SPieter Ghysels } 137d6ea485SPieter Ghysels 14291d0ed5SPieter Ghysels #undef __FUNCT__ 15291d0ed5SPieter Ghysels #define __FUNCT__ "MatDestroy_STRUMPACK" 16ad0c5e61SPieter Ghysels static PetscErrorCode MatDestroy_STRUMPACK(Mat A) 177d6ea485SPieter Ghysels { 180d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)A->spptr; 197d6ea485SPieter Ghysels PetscErrorCode ierr; 207d6ea485SPieter Ghysels PetscBool flg; 217d6ea485SPieter Ghysels 227d6ea485SPieter Ghysels PetscFunctionBegin; 237d6ea485SPieter Ghysels /* Deallocate STRUMPACK storage */ 240d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_destroy",STRUMPACK_destroy(S)); 257d6ea485SPieter Ghysels ierr = PetscFree(A->spptr);CHKERRQ(ierr); 267d6ea485SPieter Ghysels ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr); 277d6ea485SPieter Ghysels if (flg) { 287d6ea485SPieter Ghysels ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr); 297d6ea485SPieter Ghysels } else { 307d6ea485SPieter Ghysels ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr); 317d6ea485SPieter Ghysels } 32575704cbSPieter Ghysels 33575704cbSPieter Ghysels /* clear composed functions */ 34575704cbSPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr); 3534f43fa5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetReordering_C",NULL);CHKERRQ(ierr); 3634f43fa5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetColPerm_C",NULL);CHKERRQ(ierr); 3734f43fa5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSRelTol_C",NULL);CHKERRQ(ierr); 3834f43fa5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSAbsTol_C",NULL);CHKERRQ(ierr); 3934f43fa5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSMaxRank_C",NULL);CHKERRQ(ierr); 40*a36bf211SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSLeafSize_C",NULL);CHKERRQ(ierr); 41291d0ed5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSMinFrontSize_C",NULL);CHKERRQ(ierr); 42291d0ed5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSMinSepSize_C",NULL);CHKERRQ(ierr); 43575704cbSPieter Ghysels 44575704cbSPieter Ghysels PetscFunctionReturn(0); 45575704cbSPieter Ghysels } 46575704cbSPieter Ghysels 4734f43fa5SPieter Ghysels 48291d0ed5SPieter Ghysels #undef __FUNCT__ 4934f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetReordering_STRUMPACK" 5034f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetReordering_STRUMPACK(Mat F,MatSTRUMPACKReordering reordering) 51575704cbSPieter Ghysels { 520d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 53575704cbSPieter Ghysels 54575704cbSPieter Ghysels PetscFunctionBegin; 5534f43fa5SPieter Ghysels PetscStackCall("STRUMPACK_reordering_method",STRUMPACK_set_reordering_method(*S,(STRUMPACK_REORDERING_STRATEGY)reordering)); 5634f43fa5SPieter Ghysels PetscFunctionReturn(0); 5734f43fa5SPieter Ghysels } 5834f43fa5SPieter Ghysels #undef __FUNCT__ 5934f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetReordering" 6034f43fa5SPieter Ghysels /* 6134f43fa5SPieter Ghysels MatSTRUMPACKSetReordering - Set STRUMPACK fill-reducing reordering 6234f43fa5SPieter Ghysels 6334f43fa5SPieter Ghysels Input Parameters: 6434f43fa5SPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 6534f43fa5SPieter Ghysels - reordering - the code to be used to find the fill-reducing reordering 6634f43fa5SPieter Ghysels Possible values: NATURAL=0 METIS=1 PARMETIS=2 SCOTCH=3 PTSCOTCH=4 RCM=5 6734f43fa5SPieter Ghysels 6834f43fa5SPieter Ghysels Options Database: 6934f43fa5SPieter Ghysels . -mat_strumpack_reordering <METIS> - Sparsity reducing matrix reordering (choose one of) NATURAL METIS PARMETIS SCOTCH PTSCOTCH RCM (None) 7034f43fa5SPieter Ghysels 7134f43fa5SPieter Ghysels Level: beginner 7234f43fa5SPieter Ghysels 7334f43fa5SPieter Ghysels References: 7434f43fa5SPieter Ghysels . STRUMPACK manual 7534f43fa5SPieter Ghysels 7634f43fa5SPieter Ghysels .seealso: MatGetFactor() 7734f43fa5SPieter Ghysels */ 7834f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetReordering(Mat F,MatSTRUMPACKReordering reordering) 7934f43fa5SPieter Ghysels { 8034f43fa5SPieter Ghysels PetscErrorCode ierr; 8134f43fa5SPieter Ghysels 8234f43fa5SPieter Ghysels PetscFunctionBegin; 8334f43fa5SPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 8434f43fa5SPieter Ghysels PetscValidLogicalCollectiveEnum(F,reordering,2); 8534f43fa5SPieter Ghysels ierr = PetscTryMethod(F,"MatSTRUMPACKSetReordering_C",(Mat,MatSTRUMPACKReordering),(F,reordering));CHKERRQ(ierr); 86575704cbSPieter Ghysels PetscFunctionReturn(0); 87575704cbSPieter Ghysels } 88575704cbSPieter Ghysels 8934f43fa5SPieter Ghysels 90291d0ed5SPieter Ghysels #undef __FUNCT__ 9134f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetColPerm_STRUMPACK" 9234f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetColPerm_STRUMPACK(Mat F,PetscBool cperm) 9334f43fa5SPieter Ghysels { 9434f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 9534f43fa5SPieter Ghysels 9634f43fa5SPieter Ghysels PetscFunctionBegin; 9734f43fa5SPieter Ghysels PetscStackCall("STRUMPACK_set_mc64job",STRUMPACK_set_mc64job(*S,cperm ? 5 : 0)); 9834f43fa5SPieter Ghysels PetscFunctionReturn(0); 9934f43fa5SPieter Ghysels } 10034f43fa5SPieter Ghysels #undef __FUNCT__ 10134f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetColPerm" 102575704cbSPieter Ghysels /*@ 10334f43fa5SPieter Ghysels MatSTRUMPACKSetColPerm - Set whether STRUMPACK should try to permute the columns of the matrix in order to get a nonzero diagonal 104575704cbSPieter Ghysels Logically Collective on Mat 105575704cbSPieter Ghysels 106575704cbSPieter Ghysels Input Parameters: 107575704cbSPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 10834f43fa5SPieter Ghysels - cperm - whether or not to permute (internally) the columns of the matrix 109575704cbSPieter Ghysels 110575704cbSPieter Ghysels Options Database: 11134f43fa5SPieter Ghysels . -mat_strumpack_colperm <cperm> 112575704cbSPieter Ghysels 113575704cbSPieter Ghysels Level: beginner 114575704cbSPieter Ghysels 115575704cbSPieter Ghysels References: 116575704cbSPieter Ghysels . STRUMPACK manual 117575704cbSPieter Ghysels 118575704cbSPieter Ghysels .seealso: MatGetFactor() 119575704cbSPieter Ghysels @*/ 12034f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetColPerm(Mat F,PetscBool cperm) 121575704cbSPieter Ghysels { 122575704cbSPieter Ghysels PetscErrorCode ierr; 123575704cbSPieter Ghysels 124575704cbSPieter Ghysels PetscFunctionBegin; 125575704cbSPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 12634f43fa5SPieter Ghysels PetscValidLogicalCollectiveBool(F,cperm,2); 12734f43fa5SPieter Ghysels ierr = PetscTryMethod(F,"MatSTRUMPACKSetColPerm_C",(Mat,PetscBool),(F,cperm));CHKERRQ(ierr); 128575704cbSPieter Ghysels PetscFunctionReturn(0); 129575704cbSPieter Ghysels } 130575704cbSPieter Ghysels 131291d0ed5SPieter Ghysels #undef __FUNCT__ 13234f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSRelTol_STRUMPACK" 13334f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSRelTol_STRUMPACK(Mat F,PetscReal rtol) 13434f43fa5SPieter Ghysels { 13534f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 13634f43fa5SPieter Ghysels 13734f43fa5SPieter Ghysels PetscFunctionBegin; 138*a36bf211SPieter Ghysels PetscStackCall("STRUMPACK_set_HSS_rel_tol", STRUMPACK_set_HSS_rel_tol(*S,rtol)); 13934f43fa5SPieter Ghysels PetscFunctionReturn(0); 14034f43fa5SPieter Ghysels } 14134f43fa5SPieter Ghysels #undef __FUNCT__ 14234f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSRelTol" 14334f43fa5SPieter Ghysels /*@ 14434f43fa5SPieter Ghysels MatSTRUMPACKSetHSSRelTol - Set STRUMPACK relative tolerance for HSS compression 14534f43fa5SPieter Ghysels Logically Collective on Mat 14634f43fa5SPieter Ghysels 14734f43fa5SPieter Ghysels Input Parameters: 14834f43fa5SPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 14934f43fa5SPieter Ghysels - rtol - relative compression tolerance 15034f43fa5SPieter Ghysels 15134f43fa5SPieter Ghysels Options Database: 15234f43fa5SPieter Ghysels . -mat_strumpack_hss_rel_tol <1e-2> - Relative compression tolerance (None) 15334f43fa5SPieter Ghysels 15434f43fa5SPieter Ghysels Level: beginner 15534f43fa5SPieter Ghysels 15634f43fa5SPieter Ghysels References: 15734f43fa5SPieter Ghysels . STRUMPACK manual 15834f43fa5SPieter Ghysels 15934f43fa5SPieter Ghysels .seealso: MatGetFactor() 16034f43fa5SPieter Ghysels @*/ 16134f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSRelTol(Mat F,PetscReal rtol) 16234f43fa5SPieter Ghysels { 16334f43fa5SPieter Ghysels PetscErrorCode ierr; 16434f43fa5SPieter Ghysels 16534f43fa5SPieter Ghysels PetscFunctionBegin; 16634f43fa5SPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 16734f43fa5SPieter Ghysels PetscValidLogicalCollectiveReal(F,rtol,2); 16834f43fa5SPieter Ghysels ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSRelTol_C",(Mat,PetscReal),(F,rtol));CHKERRQ(ierr); 16934f43fa5SPieter Ghysels PetscFunctionReturn(0); 17034f43fa5SPieter Ghysels } 17134f43fa5SPieter Ghysels 17234f43fa5SPieter Ghysels 17334f43fa5SPieter Ghysels #undef __FUNCT__ 17434f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSAbsTol_STRUMPACK" 17534f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSAbsTol_STRUMPACK(Mat F,PetscReal atol) 17634f43fa5SPieter Ghysels { 17734f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 17834f43fa5SPieter Ghysels 17934f43fa5SPieter Ghysels PetscFunctionBegin; 180*a36bf211SPieter Ghysels PetscStackCall("STRUMPACK_set_HSS_abs_tol", STRUMPACK_set_HSS_abs_tol(*S,atol)); 18134f43fa5SPieter Ghysels PetscFunctionReturn(0); 18234f43fa5SPieter Ghysels } 18334f43fa5SPieter Ghysels #undef __FUNCT__ 18434f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSAbsTol" 18534f43fa5SPieter Ghysels /*@ 18634f43fa5SPieter Ghysels MatSTRUMPACKSetHSSAbsTol - Set STRUMPACK absolute tolerance for HSS compression 18734f43fa5SPieter Ghysels Logically Collective on Mat 18834f43fa5SPieter Ghysels 18934f43fa5SPieter Ghysels Input Parameters: 19034f43fa5SPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 19134f43fa5SPieter Ghysels - atol - absolute compression tolerance 19234f43fa5SPieter Ghysels 19334f43fa5SPieter Ghysels Options Database: 19434f43fa5SPieter Ghysels . -mat_strumpack_hss_abs_tol <1e-8> - Absolute compression tolerance (None) 19534f43fa5SPieter Ghysels 19634f43fa5SPieter Ghysels Level: beginner 19734f43fa5SPieter Ghysels 19834f43fa5SPieter Ghysels References: 19934f43fa5SPieter Ghysels . STRUMPACK manual 20034f43fa5SPieter Ghysels 20134f43fa5SPieter Ghysels .seealso: MatGetFactor() 20234f43fa5SPieter Ghysels @*/ 20334f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSAbsTol(Mat F,PetscReal atol) 20434f43fa5SPieter Ghysels { 20534f43fa5SPieter Ghysels PetscErrorCode ierr; 20634f43fa5SPieter Ghysels 20734f43fa5SPieter Ghysels PetscFunctionBegin; 20834f43fa5SPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 20934f43fa5SPieter Ghysels PetscValidLogicalCollectiveReal(F,atol,2); 21034f43fa5SPieter Ghysels ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSAbsTol_C",(Mat,PetscReal),(F,atol));CHKERRQ(ierr); 21134f43fa5SPieter Ghysels PetscFunctionReturn(0); 21234f43fa5SPieter Ghysels } 21334f43fa5SPieter Ghysels 21434f43fa5SPieter Ghysels 21534f43fa5SPieter Ghysels #undef __FUNCT__ 21634f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSMaxRank_STRUMPACK" 21734f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSMaxRank_STRUMPACK(Mat F,PetscInt hssmaxrank) 21834f43fa5SPieter Ghysels { 21934f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 22034f43fa5SPieter Ghysels 22134f43fa5SPieter Ghysels PetscFunctionBegin; 222*a36bf211SPieter Ghysels PetscStackCall("STRUMPACK_set_HSS_max_rank", STRUMPACK_set_HSS_max_rank(*S,hssmaxrank)); 22334f43fa5SPieter Ghysels PetscFunctionReturn(0); 22434f43fa5SPieter Ghysels } 22534f43fa5SPieter Ghysels #undef __FUNCT__ 22634f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSMaxRank" 22734f43fa5SPieter Ghysels /*@ 22834f43fa5SPieter Ghysels MatSTRUMPACKSetHSSMaxRank - Set STRUMPACK maximum HSS rank 22934f43fa5SPieter Ghysels Logically Collective on Mat 23034f43fa5SPieter Ghysels 23134f43fa5SPieter Ghysels Input Parameters: 23234f43fa5SPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 23334f43fa5SPieter Ghysels - hssmaxrank - maximum rank used in low-rank approximation 23434f43fa5SPieter Ghysels 23534f43fa5SPieter Ghysels Options Database: 23634f43fa5SPieter Ghysels . -mat_strumpack_max_rank - Maximum rank in HSS compression, when using pctype ilu (None) 23734f43fa5SPieter Ghysels 23834f43fa5SPieter Ghysels Level: beginner 23934f43fa5SPieter Ghysels 24034f43fa5SPieter Ghysels References: 24134f43fa5SPieter Ghysels . STRUMPACK manual 24234f43fa5SPieter Ghysels 24334f43fa5SPieter Ghysels .seealso: MatGetFactor() 24434f43fa5SPieter Ghysels @*/ 24534f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSMaxRank(Mat F,PetscInt hssmaxrank) 24634f43fa5SPieter Ghysels { 24734f43fa5SPieter Ghysels PetscErrorCode ierr; 24834f43fa5SPieter Ghysels 24934f43fa5SPieter Ghysels PetscFunctionBegin; 25034f43fa5SPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 25134f43fa5SPieter Ghysels PetscValidLogicalCollectiveInt(F,hssmaxrank,2); 25234f43fa5SPieter Ghysels ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSMaxRank_C",(Mat,PetscInt),(F,hssmaxrank));CHKERRQ(ierr); 25334f43fa5SPieter Ghysels PetscFunctionReturn(0); 25434f43fa5SPieter Ghysels } 25534f43fa5SPieter Ghysels 25634f43fa5SPieter Ghysels 25734f43fa5SPieter Ghysels #undef __FUNCT__ 258*a36bf211SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSLeafSize_STRUMPACK" 259*a36bf211SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSLeafSize_STRUMPACK(Mat F,PetscInt leaf_size) 260*a36bf211SPieter Ghysels { 261*a36bf211SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 262*a36bf211SPieter Ghysels 263*a36bf211SPieter Ghysels PetscFunctionBegin; 264*a36bf211SPieter Ghysels PetscStackCall("STRUMPACK_set_HSS_leaf_size", STRUMPACK_set_HSS_leaf_size(*S,leaf_size)); 265*a36bf211SPieter Ghysels PetscFunctionReturn(0); 266*a36bf211SPieter Ghysels } 267*a36bf211SPieter Ghysels #undef __FUNCT__ 268*a36bf211SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSLeafSize" 269*a36bf211SPieter Ghysels /*@ 270*a36bf211SPieter Ghysels MatSTRUMPACKSetHSSLeafSize - Set STRUMPACK HSS leaf size 271*a36bf211SPieter Ghysels Logically Collective on Mat 272*a36bf211SPieter Ghysels 273*a36bf211SPieter Ghysels Input Parameters: 274*a36bf211SPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 275*a36bf211SPieter Ghysels - leaf_size - Size of diagonal blocks in HSS approximation 276*a36bf211SPieter Ghysels 277*a36bf211SPieter Ghysels Options Database: 278*a36bf211SPieter Ghysels . -mat_strumpack_leaf_size - Size of diagonal blocks in HSS approximation, when using pctype ilu (None) 279*a36bf211SPieter Ghysels 280*a36bf211SPieter Ghysels Level: beginner 281*a36bf211SPieter Ghysels 282*a36bf211SPieter Ghysels References: 283*a36bf211SPieter Ghysels . STRUMPACK manual 284*a36bf211SPieter Ghysels 285*a36bf211SPieter Ghysels .seealso: MatGetFactor() 286*a36bf211SPieter Ghysels @*/ 287*a36bf211SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSLeafSize(Mat F,PetscInt leaf_size) 288*a36bf211SPieter Ghysels { 289*a36bf211SPieter Ghysels PetscErrorCode ierr; 290*a36bf211SPieter Ghysels 291*a36bf211SPieter Ghysels PetscFunctionBegin; 292*a36bf211SPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 293*a36bf211SPieter Ghysels PetscValidLogicalCollectiveInt(F,leaf_size,2); 294*a36bf211SPieter Ghysels ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSLeafSize_C",(Mat,PetscInt),(F,leaf_size));CHKERRQ(ierr); 295*a36bf211SPieter Ghysels PetscFunctionReturn(0); 296*a36bf211SPieter Ghysels } 297*a36bf211SPieter Ghysels 298*a36bf211SPieter Ghysels 299*a36bf211SPieter Ghysels #undef __FUNCT__ 300291d0ed5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSMinFrontSize_STRUMPACK" 301291d0ed5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSMinFrontSize_STRUMPACK(Mat F,PetscInt hssminsize) 302575704cbSPieter Ghysels { 3030d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 304575704cbSPieter Ghysels 305575704cbSPieter Ghysels PetscFunctionBegin; 306291d0ed5SPieter Ghysels PetscStackCall("STRUMPACK_set_HSS_min_front_size", STRUMPACK_set_HSS_min_front_size(*S,hssminsize)); 307291d0ed5SPieter Ghysels PetscFunctionReturn(0); 308291d0ed5SPieter Ghysels } 309291d0ed5SPieter Ghysels #undef __FUNCT__ 310291d0ed5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSMinFrontSize" 311575704cbSPieter Ghysels /*@ 312291d0ed5SPieter Ghysels MatSTRUMPACKSetHSSMinFrontSize - Set STRUMPACK minimum dense matrix size for low-rank approximation 313575704cbSPieter Ghysels Logically Collective on Mat 314575704cbSPieter Ghysels 315575704cbSPieter Ghysels Input Parameters: 316575704cbSPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 317575704cbSPieter Ghysels - hssminsize - minimum dense matrix size for low-rank approximation 318575704cbSPieter Ghysels 319575704cbSPieter Ghysels Options Database: 320291d0ed5SPieter Ghysels . -mat_strumpack_hss_min_front_size <hssminsize> 321575704cbSPieter Ghysels 322575704cbSPieter Ghysels Level: beginner 323575704cbSPieter Ghysels 324575704cbSPieter Ghysels References: 325575704cbSPieter Ghysels . STRUMPACK manual 326575704cbSPieter Ghysels 327575704cbSPieter Ghysels .seealso: MatGetFactor() 328575704cbSPieter Ghysels @*/ 329291d0ed5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSMinFrontSize(Mat F,PetscInt hssminsize) 330575704cbSPieter Ghysels { 331575704cbSPieter Ghysels PetscErrorCode ierr; 332575704cbSPieter Ghysels 333575704cbSPieter Ghysels PetscFunctionBegin; 334575704cbSPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 335575704cbSPieter Ghysels PetscValidLogicalCollectiveInt(F,hssminsize,2); 336291d0ed5SPieter Ghysels ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSMinFrontSize_C",(Mat,PetscInt),(F,hssminsize));CHKERRQ(ierr); 337575704cbSPieter Ghysels PetscFunctionReturn(0); 338575704cbSPieter Ghysels } 339575704cbSPieter Ghysels 340291d0ed5SPieter Ghysels #undef __FUNCT__ 34134f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSMinSepSize_STRUMPACK" 34234f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSMinSepSize_STRUMPACK(Mat F,PetscInt hssminsize) 34334f43fa5SPieter Ghysels { 34434f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 34534f43fa5SPieter Ghysels 34634f43fa5SPieter Ghysels PetscFunctionBegin; 34734f43fa5SPieter Ghysels PetscStackCall("STRUMPACK_set_HSS_min_sep_size", STRUMPACK_set_HSS_min_sep_size(*S,hssminsize)); 34834f43fa5SPieter Ghysels PetscFunctionReturn(0); 34934f43fa5SPieter Ghysels } 35034f43fa5SPieter Ghysels #undef __FUNCT__ 351291d0ed5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSMinSepSize" 352291d0ed5SPieter Ghysels /*@ 353291d0ed5SPieter Ghysels MatSTRUMPACKSetHSSMinSepSize - Set STRUMPACK minimum separator size for low-rank approximation 354291d0ed5SPieter Ghysels Logically Collective on Mat 355291d0ed5SPieter Ghysels 356291d0ed5SPieter Ghysels Input Parameters: 357291d0ed5SPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 358291d0ed5SPieter Ghysels - hssminsize - minimum dense matrix size for low-rank approximation 359291d0ed5SPieter Ghysels 360291d0ed5SPieter Ghysels Options Database: 361291d0ed5SPieter Ghysels . -mat_strumpack_hss_min_sep_size <hssminsize> 362291d0ed5SPieter Ghysels 363291d0ed5SPieter Ghysels Level: beginner 364291d0ed5SPieter Ghysels 365291d0ed5SPieter Ghysels References: 366291d0ed5SPieter Ghysels . STRUMPACK manual 367291d0ed5SPieter Ghysels 368291d0ed5SPieter Ghysels .seealso: MatGetFactor() 369291d0ed5SPieter Ghysels @*/ 370291d0ed5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat F,PetscInt hssminsize) 371291d0ed5SPieter Ghysels { 372291d0ed5SPieter Ghysels PetscErrorCode ierr; 373291d0ed5SPieter Ghysels 374291d0ed5SPieter Ghysels PetscFunctionBegin; 375291d0ed5SPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 376291d0ed5SPieter Ghysels PetscValidLogicalCollectiveInt(F,hssminsize,2); 377291d0ed5SPieter Ghysels ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSMinSepSize_C",(Mat,PetscInt),(F,hssminsize));CHKERRQ(ierr); 378291d0ed5SPieter Ghysels PetscFunctionReturn(0); 379291d0ed5SPieter Ghysels } 380291d0ed5SPieter Ghysels 381291d0ed5SPieter Ghysels 382291d0ed5SPieter Ghysels #undef __FUNCT__ 383291d0ed5SPieter Ghysels #define __FUNCT__ "MatSolve_STRUMPACK" 384ad0c5e61SPieter Ghysels static PetscErrorCode MatSolve_STRUMPACK(Mat A,Vec b_mpi,Vec x) 3857d6ea485SPieter Ghysels { 3860d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)A->spptr; 3877d6ea485SPieter Ghysels STRUMPACK_RETURN_CODE sp_err; 3887d6ea485SPieter Ghysels PetscErrorCode ierr; 3897d6ea485SPieter Ghysels const PetscScalar *bptr; 3907d6ea485SPieter Ghysels PetscScalar *xptr; 3917d6ea485SPieter Ghysels 3927d6ea485SPieter Ghysels PetscFunctionBegin; 3937d6ea485SPieter Ghysels ierr = VecGetArray(x,&xptr);CHKERRQ(ierr); 3947d6ea485SPieter Ghysels ierr = VecGetArrayRead(b_mpi,&bptr);CHKERRQ(ierr); 3950d08a34dSPieter Ghysels 3960d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_solve",sp_err = STRUMPACK_solve(*S,(PetscScalar*)bptr,xptr,0)); 3970d08a34dSPieter Ghysels switch (sp_err) { 3980d08a34dSPieter Ghysels case STRUMPACK_SUCCESS: break; 3990d08a34dSPieter Ghysels case STRUMPACK_MATRIX_NOT_SET: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); break; } 4000d08a34dSPieter Ghysels case STRUMPACK_REORDERING_ERROR: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); break; } 4010d08a34dSPieter Ghysels default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: solve failed"); 4027d6ea485SPieter Ghysels } 4037d6ea485SPieter Ghysels ierr = VecRestoreArray(x,&xptr);CHKERRQ(ierr); 4047d6ea485SPieter Ghysels ierr = VecRestoreArrayRead(b_mpi,&bptr);CHKERRQ(ierr); 4057d6ea485SPieter Ghysels PetscFunctionReturn(0); 4067d6ea485SPieter Ghysels } 4077d6ea485SPieter Ghysels 408291d0ed5SPieter Ghysels #undef __FUNCT__ 409291d0ed5SPieter Ghysels #define __FUNCT__ "MatMatSolve_STRUMPACK" 410ad0c5e61SPieter Ghysels static PetscErrorCode MatMatSolve_STRUMPACK(Mat A,Mat B_mpi,Mat X) 4117d6ea485SPieter Ghysels { 4127d6ea485SPieter Ghysels PetscErrorCode ierr; 4137d6ea485SPieter Ghysels PetscBool flg; 4147d6ea485SPieter Ghysels 4157d6ea485SPieter Ghysels PetscFunctionBegin; 4167d6ea485SPieter Ghysels ierr = PetscObjectTypeCompareAny((PetscObject)B_mpi,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 4177d6ea485SPieter Ghysels if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 4187d6ea485SPieter Ghysels ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 4197d6ea485SPieter Ghysels if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); 4207d6ea485SPieter Ghysels SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_STRUMPACK() is not implemented yet"); 4217d6ea485SPieter Ghysels PetscFunctionReturn(0); 4227d6ea485SPieter Ghysels } 4237d6ea485SPieter Ghysels 424291d0ed5SPieter Ghysels #undef __FUNCT__ 425291d0ed5SPieter Ghysels #define __FUNCT__ "MatFactorInfo_STRUMPACK" 426ad0c5e61SPieter Ghysels static PetscErrorCode MatFactorInfo_STRUMPACK(Mat A,PetscViewer viewer) 427ad0c5e61SPieter Ghysels { 428ad0c5e61SPieter Ghysels PetscErrorCode ierr; 429ad0c5e61SPieter Ghysels 430ad0c5e61SPieter Ghysels PetscFunctionBegin; 431ad0c5e61SPieter Ghysels /* check if matrix is strumpack type */ 432ad0c5e61SPieter Ghysels if (A->ops->solve != MatSolve_STRUMPACK) PetscFunctionReturn(0); 433ad0c5e61SPieter Ghysels ierr = PetscViewerASCIIPrintf(viewer,"STRUMPACK sparse solver!\n");CHKERRQ(ierr); 434ad0c5e61SPieter Ghysels PetscFunctionReturn(0); 435ad0c5e61SPieter Ghysels } 436ad0c5e61SPieter Ghysels 437291d0ed5SPieter Ghysels #undef __FUNCT__ 438291d0ed5SPieter Ghysels #define __FUNCT__ "MatView_STRUMPACK" 439ad0c5e61SPieter Ghysels static PetscErrorCode MatView_STRUMPACK(Mat A,PetscViewer viewer) 440ad0c5e61SPieter Ghysels { 441ad0c5e61SPieter Ghysels PetscErrorCode ierr; 442ad0c5e61SPieter Ghysels PetscBool iascii; 443ad0c5e61SPieter Ghysels PetscViewerFormat format; 444ad0c5e61SPieter Ghysels 445ad0c5e61SPieter Ghysels PetscFunctionBegin; 446ad0c5e61SPieter Ghysels ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 447ad0c5e61SPieter Ghysels if (iascii) { 448ad0c5e61SPieter Ghysels ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 449ad0c5e61SPieter Ghysels if (format == PETSC_VIEWER_ASCII_INFO) { 450ad0c5e61SPieter Ghysels ierr = MatFactorInfo_STRUMPACK(A,viewer);CHKERRQ(ierr); 451ad0c5e61SPieter Ghysels } 452ad0c5e61SPieter Ghysels } 453ad0c5e61SPieter Ghysels PetscFunctionReturn(0); 454ad0c5e61SPieter Ghysels } 4557d6ea485SPieter Ghysels 456291d0ed5SPieter Ghysels #undef __FUNCT__ 457291d0ed5SPieter Ghysels #define __FUNCT__ "MatLUFactorNumeric_STRUMPACK" 458ad0c5e61SPieter Ghysels static PetscErrorCode MatLUFactorNumeric_STRUMPACK(Mat F,Mat A,const MatFactorInfo *info) 4597d6ea485SPieter Ghysels { 4600d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 4617d6ea485SPieter Ghysels STRUMPACK_RETURN_CODE sp_err; 4627d6ea485SPieter Ghysels Mat_SeqAIJ *A_d,*A_o; 4637d6ea485SPieter Ghysels Mat_MPIAIJ *mat; 4647d6ea485SPieter Ghysels PetscErrorCode ierr; 4650d08a34dSPieter Ghysels PetscInt M=A->rmap->N,m=A->rmap->n; 4667d6ea485SPieter Ghysels PetscBool flg; 4677d6ea485SPieter Ghysels 4687d6ea485SPieter Ghysels PetscFunctionBegin; 4697d6ea485SPieter Ghysels ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); 4707d6ea485SPieter Ghysels if (flg) { /* A is MATMPIAIJ */ 4717d6ea485SPieter Ghysels mat = (Mat_MPIAIJ*)A->data; 4727d6ea485SPieter Ghysels A_d = (Mat_SeqAIJ*)(mat->A)->data; 4737d6ea485SPieter Ghysels A_o = (Mat_SeqAIJ*)(mat->B)->data; 4740d08a34dSPieter Ghysels PetscStackCall("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)); 4757d6ea485SPieter Ghysels } else { /* A is MATSEQAIJ */ 4767d6ea485SPieter Ghysels A_d = (Mat_SeqAIJ*)A->data; 4770d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_set_csr_matrix",STRUMPACK_set_csr_matrix(*S,&M,A_d->i,A_d->j,A_d->a,0)); 4787d6ea485SPieter Ghysels } 4797d6ea485SPieter Ghysels 4807d6ea485SPieter Ghysels /* Reorder and Factor the matrix. */ 4817d6ea485SPieter Ghysels /* TODO figure out how to avoid reorder if the matrix values changed, but the pattern remains the same. */ 4820d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_reorder",sp_err = STRUMPACK_reorder(*S)); 4830d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_factor",sp_err = STRUMPACK_factor(*S)); 4840d08a34dSPieter Ghysels switch (sp_err) { 4850d08a34dSPieter Ghysels case STRUMPACK_SUCCESS: break; 4860d08a34dSPieter Ghysels case STRUMPACK_MATRIX_NOT_SET: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); break; } 4870d08a34dSPieter Ghysels case STRUMPACK_REORDERING_ERROR: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); break; } 4880d08a34dSPieter Ghysels default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: factorization failed"); 4897d6ea485SPieter Ghysels } 4907d6ea485SPieter Ghysels PetscFunctionReturn(0); 4917d6ea485SPieter Ghysels } 4927d6ea485SPieter Ghysels 493291d0ed5SPieter Ghysels #undef __FUNCT__ 494291d0ed5SPieter Ghysels #define __FUNCT__ "MatLUFactorSymbolic_STRUMPACK" 495ad0c5e61SPieter Ghysels static PetscErrorCode MatLUFactorSymbolic_STRUMPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 4967d6ea485SPieter Ghysels { 4977d6ea485SPieter Ghysels PetscFunctionBegin; 4987d6ea485SPieter Ghysels F->ops->lufactornumeric = MatLUFactorNumeric_STRUMPACK; 4997d6ea485SPieter Ghysels F->ops->solve = MatSolve_STRUMPACK; 5007d6ea485SPieter Ghysels F->ops->matsolve = MatMatSolve_STRUMPACK; 5017d6ea485SPieter Ghysels PetscFunctionReturn(0); 5027d6ea485SPieter Ghysels } 5037d6ea485SPieter Ghysels 504291d0ed5SPieter Ghysels #undef __FUNCT__ 505291d0ed5SPieter Ghysels #define __FUNCT__ "MatFactorGetSolverPackage_aij_strumpack" 506ad0c5e61SPieter Ghysels static PetscErrorCode MatFactorGetSolverPackage_aij_strumpack(Mat A,const MatSolverPackage *type) 5077d6ea485SPieter Ghysels { 5087d6ea485SPieter Ghysels PetscFunctionBegin; 5097d6ea485SPieter Ghysels *type = MATSOLVERSTRUMPACK; 5107d6ea485SPieter Ghysels PetscFunctionReturn(0); 5117d6ea485SPieter Ghysels } 5127d6ea485SPieter Ghysels 513575704cbSPieter Ghysels /*MC 514575704cbSPieter Ghysels MATSOLVERSSTRUMPACK = "strumpack" - A solver package providing a direct sparse solver (PCLU) 515575704cbSPieter Ghysels and a preconditioner (PCILU) using low-rank compression via the external package STRUMPACK. 516575704cbSPieter Ghysels 517575704cbSPieter Ghysels Consult the STRUMPACK-sparse manual for more info. 518575704cbSPieter Ghysels 519575704cbSPieter Ghysels Use 520575704cbSPieter Ghysels ./configure --download-strumpack 521575704cbSPieter Ghysels to have PETSc installed with STRUMPACK 522575704cbSPieter Ghysels 523575704cbSPieter Ghysels Use 524575704cbSPieter Ghysels -pc_type lu -pc_factor_mat_solver_package strumpack 525575704cbSPieter Ghysels to use this as an exact (direct) solver, use 526575704cbSPieter Ghysels -pc_type ilu -pc_factor_mat_solver_package strumpack 527575704cbSPieter Ghysels to enable low-rank compression (i.e, use as a preconditioner). 528575704cbSPieter Ghysels 529575704cbSPieter Ghysels Works with AIJ matrices 530575704cbSPieter Ghysels 531575704cbSPieter Ghysels Options Database Keys: 53234f43fa5SPieter Ghysels + -mat_strumpack_verbose 53334f43fa5SPieter Ghysels . -mat_strumpack_hss_rel_tol <1e-2> - Relative compression tolerance (None) 53434f43fa5SPieter Ghysels . -mat_strumpack_hss_abs_tol <1e-8> - Absolute compression tolerance (None) 535575704cbSPieter Ghysels . -mat_strumpack_colperm <TRUE> - Permute matrix to make diagonal nonzeros (None) 53634f43fa5SPieter Ghysels . -mat_strumpack_hss_min_front_size <1000> - Minimum size of dense block for HSS compression (None) 53734f43fa5SPieter Ghysels . -mat_strumpack_hss_min_sep_size <256> - Minimum size of separator for HSS compression (None) 53834f43fa5SPieter Ghysels . -mat_strumpack_max_rank - Maximum rank in HSS compression, when using pctype ilu (None) 539*a36bf211SPieter Ghysels . -mat_strumpack_leaf_size - Size of diagonal blocks in HSS approximation, when using pctype ilu (None) 54034f43fa5SPieter Ghysels . -mat_strumpack_reordering <METIS> - Sparsity reducing matrix reordering (choose one of) NATURAL METIS PARMETIS SCOTCH PTSCOTCH RCM (None) 54134f43fa5SPieter Ghysels . -mat_strumpack_iterative_solver <DIRECT> - Select iterative solver from STRUMPACK (choose one of) AUTO DIRECT REFINE PREC_GMRES GMRES PREC_BICGSTAB BICGSTAB (None) 542575704cbSPieter Ghysels 543575704cbSPieter Ghysels Level: beginner 544575704cbSPieter Ghysels 545575704cbSPieter Ghysels .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, PCFactorSetMatSolverPackage(), MatSolverPackage 546575704cbSPieter Ghysels M*/ 547291d0ed5SPieter Ghysels #undef __FUNCT__ 548291d0ed5SPieter Ghysels #define __FUNCT__ "MatGetFactor_aij_strumpack" 549ad0c5e61SPieter Ghysels static PetscErrorCode MatGetFactor_aij_strumpack(Mat A,MatFactorType ftype,Mat *F) 5507d6ea485SPieter Ghysels { 5517d6ea485SPieter Ghysels Mat B; 5527d6ea485SPieter Ghysels PetscErrorCode ierr; 5537d6ea485SPieter Ghysels PetscInt M=A->rmap->N,N=A->cmap->N; 554575704cbSPieter Ghysels PetscBool verb,flg,set; 55534f43fa5SPieter Ghysels PetscReal ctol; 556*a36bf211SPieter Ghysels PetscInt hssminsize,max_rank,leaf_size; 55734f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S; 55834f43fa5SPieter Ghysels STRUMPACK_INTERFACE iface; 55934f43fa5SPieter Ghysels STRUMPACK_REORDERING_STRATEGY ndcurrent,ndvalue; 56034f43fa5SPieter Ghysels STRUMPACK_KRYLOV_SOLVER itcurrent,itsolver; 56134f43fa5SPieter Ghysels const STRUMPACK_PRECISION table[2][2][2] = 56234f43fa5SPieter Ghysels {{{STRUMPACK_FLOATCOMPLEX_64, STRUMPACK_DOUBLECOMPLEX_64}, 5637d6ea485SPieter Ghysels {STRUMPACK_FLOAT_64, STRUMPACK_DOUBLE_64}}, 5647d6ea485SPieter Ghysels {{STRUMPACK_FLOATCOMPLEX, STRUMPACK_DOUBLECOMPLEX}, 5657d6ea485SPieter Ghysels {STRUMPACK_FLOAT, STRUMPACK_DOUBLE}}}; 56634f43fa5SPieter Ghysels const STRUMPACK_PRECISION prec = 56734f43fa5SPieter Ghysels table[(sizeof(PetscInt)==8)?0:1] 56834f43fa5SPieter Ghysels [(PETSC_SCALAR==PETSC_COMPLEX)?0:1] 56934f43fa5SPieter Ghysels [(PETSC_REAL==PETSC_FLOAT)?0:1]; 5707d6ea485SPieter Ghysels 5717d6ea485SPieter Ghysels PetscFunctionBegin; 5727d6ea485SPieter Ghysels /* Create the factorization matrix */ 5737d6ea485SPieter Ghysels ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 5747d6ea485SPieter Ghysels ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,M,N);CHKERRQ(ierr); 5757d6ea485SPieter Ghysels ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 5767d6ea485SPieter Ghysels ierr = MatSeqAIJSetPreallocation(B,0,NULL); 5777d6ea485SPieter Ghysels ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr); 578575704cbSPieter Ghysels if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) { 5797d6ea485SPieter Ghysels B->ops->lufactorsymbolic = MatLUFactorSymbolic_STRUMPACK; 580575704cbSPieter Ghysels B->ops->ilufactorsymbolic = MatLUFactorSymbolic_STRUMPACK; 581575704cbSPieter Ghysels } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 5827d6ea485SPieter Ghysels B->ops->view = MatView_STRUMPACK; 5837d6ea485SPieter Ghysels B->ops->destroy = MatDestroy_STRUMPACK; 5847d6ea485SPieter Ghysels B->ops->getdiagonal = MatGetDiagonal_STRUMPACK; 5857d6ea485SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_aij_strumpack);CHKERRQ(ierr); 58634f43fa5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetReordering_C",MatSTRUMPACKSetReordering_STRUMPACK);CHKERRQ(ierr); 58734f43fa5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetColPerm_C",MatSTRUMPACKSetColPerm_STRUMPACK);CHKERRQ(ierr); 58834f43fa5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSRelTol_C",MatSTRUMPACKSetHSSRelTol_STRUMPACK);CHKERRQ(ierr); 58934f43fa5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSAbsTol_C",MatSTRUMPACKSetHSSAbsTol_STRUMPACK);CHKERRQ(ierr); 590*a36bf211SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSMaxRank_C",MatSTRUMPACKSetHSSMaxRank_STRUMPACK);CHKERRQ(ierr); 591*a36bf211SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSLeafSize_C",MatSTRUMPACKSetHSSLeafSize_STRUMPACK);CHKERRQ(ierr); 592291d0ed5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSMinFrontSize_C",MatSTRUMPACKSetHSSMinFrontSize_STRUMPACK);CHKERRQ(ierr); 593291d0ed5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSMinSepSize_C",MatSTRUMPACKSetHSSMinSepSize_STRUMPACK);CHKERRQ(ierr); 594575704cbSPieter Ghysels B->factortype = ftype; 5950d08a34dSPieter Ghysels ierr = PetscNewLog(B,&S);CHKERRQ(ierr); 5960d08a34dSPieter Ghysels B->spptr = S; 5970d08a34dSPieter Ghysels 5980d08a34dSPieter Ghysels ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg); 5990d08a34dSPieter Ghysels iface = flg ? STRUMPACK_MT : STRUMPACK_MPI_DIST; 6007d6ea485SPieter Ghysels 6017d6ea485SPieter Ghysels ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"STRUMPACK Options","Mat");CHKERRQ(ierr); 6027d6ea485SPieter Ghysels 603575704cbSPieter Ghysels verb = PetscLogPrintInfo ? PETSC_TRUE : PETSC_FALSE; 6047d6ea485SPieter Ghysels ierr = PetscOptionsBool("-mat_strumpack_verbose","Print STRUMPACK information","None",verb,&verb,NULL);CHKERRQ(ierr); 6057d6ea485SPieter Ghysels 60634f43fa5SPieter Ghysels PetscStackCall("STRUMPACK_init",STRUMPACK_init(S,PetscObjectComm((PetscObject)A),prec,iface,0,NULL,verb)); 60755c022e5SPieter Ghysels 608*a36bf211SPieter Ghysels PetscStackCall("STRUMPACK_HSS_rel_tol",ctol = (PetscReal)STRUMPACK_HSS_rel_tol(*S)); 60934f43fa5SPieter Ghysels ierr = PetscOptionsReal("-mat_strumpack_hss_rel_tol","Relative compression tolerance","None",ctol,&ctol,&set);CHKERRQ(ierr); 610*a36bf211SPieter Ghysels if (set) PetscStackCall("STRUMPACK_set_HSS_rel_tol",STRUMPACK_set_HSS_rel_tol(*S,(double)ctol)); 6117d6ea485SPieter Ghysels 612*a36bf211SPieter Ghysels PetscStackCall("STRUMPACK_HSS_abs_tol",ctol = (PetscReal)STRUMPACK_HSS_abs_tol(*S)); 61334f43fa5SPieter Ghysels ierr = PetscOptionsReal("-mat_strumpack_hss_abs_tol","Absolute compression tolerance","None",ctol,&ctol,&set);CHKERRQ(ierr); 614*a36bf211SPieter Ghysels if (set) PetscStackCall("STRUMPACK_set_HSS_abs_tol",STRUMPACK_set_HSS_abs_tol(*S,(double)ctol)); 615575704cbSPieter Ghysels 616291d0ed5SPieter Ghysels PetscStackCall("STRUMPACK_mc64job",flg = (STRUMPACK_mc64job(*S) == 0) ? PETSC_FALSE : PETSC_TRUE); 617575704cbSPieter Ghysels ierr = PetscOptionsBool("-mat_strumpack_colperm","Find a col perm to get nonzero diagonal","None",flg,&flg,&set);CHKERRQ(ierr); 6180d08a34dSPieter Ghysels if (set) PetscStackCall("STRUMPACK_set_mc64job",STRUMPACK_set_mc64job(*S,flg ? 5 : 0)); 619575704cbSPieter Ghysels 620291d0ed5SPieter Ghysels PetscStackCall("STRUMPACK_HSS_min_front_size",hssminsize = (PetscInt)STRUMPACK_HSS_min_front_size(*S)); 621291d0ed5SPieter Ghysels ierr = PetscOptionsInt("-mat_strumpack_hss_min_front_size","Minimum size of dense block for HSS compression","None",hssminsize,&hssminsize,&set);CHKERRQ(ierr); 622*a36bf211SPieter Ghysels if (set) { 623*a36bf211SPieter Ghysels printf("setting HSS min front size to %d\n", hssminsize); 624*a36bf211SPieter Ghysels PetscStackCall("STRUMPACK_set_HSS_min_front_size",STRUMPACK_set_HSS_min_front_size(*S,(int)hssminsize)); 625*a36bf211SPieter Ghysels } 626291d0ed5SPieter Ghysels 627291d0ed5SPieter Ghysels PetscStackCall("STRUMPACK_HSS_min_sep_size",hssminsize = (PetscInt)STRUMPACK_HSS_min_sep_size(*S)); 628291d0ed5SPieter Ghysels ierr = PetscOptionsInt("-mat_strumpack_hss_min_sep_size","Minimum size of separator for HSS compression","None",hssminsize,&hssminsize,&set);CHKERRQ(ierr); 629291d0ed5SPieter Ghysels if (set) PetscStackCall("STRUMPACK_set_HSS_min_sep_size",STRUMPACK_set_HSS_min_sep_size(*S,(int)hssminsize)); 630575704cbSPieter Ghysels 631*a36bf211SPieter Ghysels PetscStackCall("STRUMPACK_HSS_max_rank",max_rank = (PetscInt)STRUMPACK_HSS_max_rank(*S)); 63234f43fa5SPieter Ghysels ierr = PetscOptionsInt("-mat_strumpack_max_rank","Maximum rank in HSS approximation","None",max_rank,&max_rank,&set);CHKERRQ(ierr); 633*a36bf211SPieter Ghysels if (set) PetscStackCall("STRUMPACK_set_HSS_max_rank",STRUMPACK_set_HSS_max_rank(*S,(int)max_rank)); 634*a36bf211SPieter Ghysels 635*a36bf211SPieter Ghysels PetscStackCall("STRUMPACK_HSS_leaf_size",leaf_size = (PetscInt)STRUMPACK_HSS_leaf_size(*S)); 636*a36bf211SPieter Ghysels ierr = PetscOptionsInt("-mat_strumpack_leaf_size","Size of diagonal blocks in HSS approximation","None",leaf_size,&leaf_size,&set);CHKERRQ(ierr); 637*a36bf211SPieter Ghysels if (set) PetscStackCall("STRUMPACK_set_HSS_leaf_size",STRUMPACK_set_HSS_leaf_size(*S,(int)leaf_size)); 63834f43fa5SPieter Ghysels 63934f43fa5SPieter Ghysels const char *const STRUMPACKNDTypes[] = {"NATURAL","METIS","PARMETIS","SCOTCH","PTSCOTCH","RCM","STRUMPACKNDTypes","",0}; 64034f43fa5SPieter Ghysels PetscStackCall("STRUMPACK_reordering_method",ndcurrent = STRUMPACK_reordering_method(*S)); 64134f43fa5SPieter Ghysels PetscOptionsEnum("-mat_strumpack_reordering","Sparsity reducing matrix reordering","None",STRUMPACKNDTypes,(PetscEnum)ndcurrent,(PetscEnum*)&ndvalue,&set);CHKERRQ(ierr); 64234f43fa5SPieter Ghysels if (set) PetscStackCall("STRUMPACK_set_reordering_method",STRUMPACK_set_reordering_method(*S,ndvalue)); 643575704cbSPieter Ghysels 64488382b05SPieter Ghysels if (ftype == MAT_FACTOR_ILU) { 64588382b05SPieter Ghysels /* When enabling HSS compression, the STRUMPACK solver becomes an incomplete */ 64688382b05SPieter Ghysels /* (or approximate) LU factorization. */ 647291d0ed5SPieter Ghysels PetscStackCall("STRUMPACK_enable_HSS",STRUMPACK_enable_HSS(*S)); 64888382b05SPieter Ghysels } 649575704cbSPieter Ghysels 65034f43fa5SPieter Ghysels /* Disable the outer iterative solver from STRUMPACK. */ 65134f43fa5SPieter Ghysels /* When STRUMPACK is used as a direct solver, it will by default do iterative refinement. */ 65234f43fa5SPieter Ghysels /* When STRUMPACK is used as an approximate factorization preconditioner (by enabling */ 65334f43fa5SPieter Ghysels /* low-rank compression), it will use it's own preconditioned GMRES. Here we can disable */ 65434f43fa5SPieter Ghysels /* the outer iterative solver, as PETSc uses STRUMPACK from within a KSP. */ 65534f43fa5SPieter Ghysels PetscStackCall("STRUMPACK_set_Krylov_solver", STRUMPACK_set_Krylov_solver(*S, STRUMPACK_DIRECT)); 65634f43fa5SPieter Ghysels 65734f43fa5SPieter Ghysels const char *const SolverTypes[] = {"AUTO","NONE","REFINE","PREC_GMRES","GMRES","PREC_BICGSTAB","BICGSTAB","SolverTypes","",0}; 65834f43fa5SPieter Ghysels PetscStackCall("STRUMPACK_Krylov_solver",itcurrent = STRUMPACK_Krylov_solver(*S)); 65934f43fa5SPieter Ghysels PetscOptionsEnum("-mat_strumpack_iterative_solver","Select iterative solver from STRUMPACK","None",SolverTypes,(PetscEnum)itcurrent,(PetscEnum*)&itsolver,&set);CHKERRQ(ierr); 66034f43fa5SPieter Ghysels if (set) PetscStackCall("STRUMPACK_set_Krylov_solver",STRUMPACK_set_Krylov_solver(*S,itsolver)); 66134f43fa5SPieter Ghysels 66234f43fa5SPieter Ghysels PetscOptionsEnd(); 66355c022e5SPieter Ghysels 6647d6ea485SPieter Ghysels *F = B; 6657d6ea485SPieter Ghysels PetscFunctionReturn(0); 6667d6ea485SPieter Ghysels } 6677d6ea485SPieter Ghysels 668291d0ed5SPieter Ghysels #undef __FUNCT__ 669291d0ed5SPieter Ghysels #define __FUNCT__ "MatSolverPackageRegister_STRUMPACK" 6707d6ea485SPieter Ghysels PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_STRUMPACK(void) 6717d6ea485SPieter Ghysels { 6727d6ea485SPieter Ghysels PetscErrorCode ierr; 6737d6ea485SPieter Ghysels 6747d6ea485SPieter Ghysels PetscFunctionBegin; 6757d6ea485SPieter Ghysels ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr); 6767d6ea485SPieter Ghysels ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr); 677575704cbSPieter Ghysels ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_ILU,MatGetFactor_aij_strumpack);CHKERRQ(ierr); 678575704cbSPieter Ghysels ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_ILU,MatGetFactor_aij_strumpack);CHKERRQ(ierr); 6797d6ea485SPieter Ghysels PetscFunctionReturn(0); 6807d6ea485SPieter Ghysels } 681