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); 35*34f43fa5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetReordering_C",NULL);CHKERRQ(ierr); 36*34f43fa5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetColPerm_C",NULL);CHKERRQ(ierr); 37*34f43fa5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSRelTol_C",NULL);CHKERRQ(ierr); 38*34f43fa5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSAbsTol_C",NULL);CHKERRQ(ierr); 39*34f43fa5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSMaxRank_C",NULL);CHKERRQ(ierr); 40291d0ed5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSMinFrontSize_C",NULL);CHKERRQ(ierr); 41291d0ed5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSMinSepSize_C",NULL);CHKERRQ(ierr); 42575704cbSPieter Ghysels 43575704cbSPieter Ghysels PetscFunctionReturn(0); 44575704cbSPieter Ghysels } 45575704cbSPieter Ghysels 46*34f43fa5SPieter Ghysels 47291d0ed5SPieter Ghysels #undef __FUNCT__ 48*34f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetReordering_STRUMPACK" 49*34f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetReordering_STRUMPACK(Mat F,MatSTRUMPACKReordering reordering) 50575704cbSPieter Ghysels { 510d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 52575704cbSPieter Ghysels 53575704cbSPieter Ghysels PetscFunctionBegin; 54*34f43fa5SPieter Ghysels PetscStackCall("STRUMPACK_reordering_method",STRUMPACK_set_reordering_method(*S,(STRUMPACK_REORDERING_STRATEGY)reordering)); 55*34f43fa5SPieter Ghysels PetscFunctionReturn(0); 56*34f43fa5SPieter Ghysels } 57*34f43fa5SPieter Ghysels #undef __FUNCT__ 58*34f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetReordering" 59*34f43fa5SPieter Ghysels /* 60*34f43fa5SPieter Ghysels MatSTRUMPACKSetReordering - Set STRUMPACK fill-reducing reordering 61*34f43fa5SPieter Ghysels 62*34f43fa5SPieter Ghysels Input Parameters: 63*34f43fa5SPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 64*34f43fa5SPieter Ghysels - reordering - the code to be used to find the fill-reducing reordering 65*34f43fa5SPieter Ghysels Possible values: NATURAL=0 METIS=1 PARMETIS=2 SCOTCH=3 PTSCOTCH=4 RCM=5 66*34f43fa5SPieter Ghysels 67*34f43fa5SPieter Ghysels Options Database: 68*34f43fa5SPieter Ghysels . -mat_strumpack_reordering <METIS> - Sparsity reducing matrix reordering (choose one of) NATURAL METIS PARMETIS SCOTCH PTSCOTCH RCM (None) 69*34f43fa5SPieter Ghysels 70*34f43fa5SPieter Ghysels Level: beginner 71*34f43fa5SPieter Ghysels 72*34f43fa5SPieter Ghysels References: 73*34f43fa5SPieter Ghysels . STRUMPACK manual 74*34f43fa5SPieter Ghysels 75*34f43fa5SPieter Ghysels .seealso: MatGetFactor() 76*34f43fa5SPieter Ghysels */ 77*34f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetReordering(Mat F,MatSTRUMPACKReordering reordering) 78*34f43fa5SPieter Ghysels { 79*34f43fa5SPieter Ghysels PetscErrorCode ierr; 80*34f43fa5SPieter Ghysels 81*34f43fa5SPieter Ghysels PetscFunctionBegin; 82*34f43fa5SPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 83*34f43fa5SPieter Ghysels PetscValidLogicalCollectiveEnum(F,reordering,2); 84*34f43fa5SPieter Ghysels ierr = PetscTryMethod(F,"MatSTRUMPACKSetReordering_C",(Mat,MatSTRUMPACKReordering),(F,reordering));CHKERRQ(ierr); 85575704cbSPieter Ghysels PetscFunctionReturn(0); 86575704cbSPieter Ghysels } 87575704cbSPieter Ghysels 88*34f43fa5SPieter Ghysels 89291d0ed5SPieter Ghysels #undef __FUNCT__ 90*34f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetColPerm_STRUMPACK" 91*34f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetColPerm_STRUMPACK(Mat F,PetscBool cperm) 92*34f43fa5SPieter Ghysels { 93*34f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 94*34f43fa5SPieter Ghysels 95*34f43fa5SPieter Ghysels PetscFunctionBegin; 96*34f43fa5SPieter Ghysels PetscStackCall("STRUMPACK_set_mc64job",STRUMPACK_set_mc64job(*S,cperm ? 5 : 0)); 97*34f43fa5SPieter Ghysels PetscFunctionReturn(0); 98*34f43fa5SPieter Ghysels } 99*34f43fa5SPieter Ghysels #undef __FUNCT__ 100*34f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetColPerm" 101575704cbSPieter Ghysels /*@ 102*34f43fa5SPieter Ghysels MatSTRUMPACKSetColPerm - Set whether STRUMPACK should try to permute the columns of the matrix in order to get a nonzero diagonal 103575704cbSPieter Ghysels Logically Collective on Mat 104575704cbSPieter Ghysels 105575704cbSPieter Ghysels Input Parameters: 106575704cbSPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 107*34f43fa5SPieter Ghysels - cperm - whether or not to permute (internally) the columns of the matrix 108575704cbSPieter Ghysels 109575704cbSPieter Ghysels Options Database: 110*34f43fa5SPieter Ghysels . -mat_strumpack_colperm <cperm> 111575704cbSPieter Ghysels 112575704cbSPieter Ghysels Level: beginner 113575704cbSPieter Ghysels 114575704cbSPieter Ghysels References: 115575704cbSPieter Ghysels . STRUMPACK manual 116575704cbSPieter Ghysels 117575704cbSPieter Ghysels .seealso: MatGetFactor() 118575704cbSPieter Ghysels @*/ 119*34f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetColPerm(Mat F,PetscBool cperm) 120575704cbSPieter Ghysels { 121575704cbSPieter Ghysels PetscErrorCode ierr; 122575704cbSPieter Ghysels 123575704cbSPieter Ghysels PetscFunctionBegin; 124575704cbSPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 125*34f43fa5SPieter Ghysels PetscValidLogicalCollectiveBool(F,cperm,2); 126*34f43fa5SPieter Ghysels ierr = PetscTryMethod(F,"MatSTRUMPACKSetColPerm_C",(Mat,PetscBool),(F,cperm));CHKERRQ(ierr); 127575704cbSPieter Ghysels PetscFunctionReturn(0); 128575704cbSPieter Ghysels } 129575704cbSPieter Ghysels 130291d0ed5SPieter Ghysels #undef __FUNCT__ 131*34f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSRelTol_STRUMPACK" 132*34f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSRelTol_STRUMPACK(Mat F,PetscReal rtol) 133*34f43fa5SPieter Ghysels { 134*34f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 135*34f43fa5SPieter Ghysels 136*34f43fa5SPieter Ghysels PetscFunctionBegin; 137*34f43fa5SPieter Ghysels PetscStackCall("STRUMPACK_set_hss_rel_tol", STRUMPACK_set_hss_rel_tol(*S,rtol)); 138*34f43fa5SPieter Ghysels PetscFunctionReturn(0); 139*34f43fa5SPieter Ghysels } 140*34f43fa5SPieter Ghysels #undef __FUNCT__ 141*34f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSRelTol" 142*34f43fa5SPieter Ghysels /*@ 143*34f43fa5SPieter Ghysels MatSTRUMPACKSetHSSRelTol - Set STRUMPACK relative tolerance for HSS compression 144*34f43fa5SPieter Ghysels Logically Collective on Mat 145*34f43fa5SPieter Ghysels 146*34f43fa5SPieter Ghysels Input Parameters: 147*34f43fa5SPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 148*34f43fa5SPieter Ghysels - rtol - relative compression tolerance 149*34f43fa5SPieter Ghysels 150*34f43fa5SPieter Ghysels Options Database: 151*34f43fa5SPieter Ghysels . -mat_strumpack_hss_rel_tol <1e-2> - Relative compression tolerance (None) 152*34f43fa5SPieter Ghysels 153*34f43fa5SPieter Ghysels Level: beginner 154*34f43fa5SPieter Ghysels 155*34f43fa5SPieter Ghysels References: 156*34f43fa5SPieter Ghysels . STRUMPACK manual 157*34f43fa5SPieter Ghysels 158*34f43fa5SPieter Ghysels .seealso: MatGetFactor() 159*34f43fa5SPieter Ghysels @*/ 160*34f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSRelTol(Mat F,PetscReal rtol) 161*34f43fa5SPieter Ghysels { 162*34f43fa5SPieter Ghysels PetscErrorCode ierr; 163*34f43fa5SPieter Ghysels 164*34f43fa5SPieter Ghysels PetscFunctionBegin; 165*34f43fa5SPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 166*34f43fa5SPieter Ghysels PetscValidLogicalCollectiveReal(F,rtol,2); 167*34f43fa5SPieter Ghysels ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSRelTol_C",(Mat,PetscReal),(F,rtol));CHKERRQ(ierr); 168*34f43fa5SPieter Ghysels PetscFunctionReturn(0); 169*34f43fa5SPieter Ghysels } 170*34f43fa5SPieter Ghysels 171*34f43fa5SPieter Ghysels 172*34f43fa5SPieter Ghysels #undef __FUNCT__ 173*34f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSAbsTol_STRUMPACK" 174*34f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSAbsTol_STRUMPACK(Mat F,PetscReal atol) 175*34f43fa5SPieter Ghysels { 176*34f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 177*34f43fa5SPieter Ghysels 178*34f43fa5SPieter Ghysels PetscFunctionBegin; 179*34f43fa5SPieter Ghysels PetscStackCall("STRUMPACK_set_hss_abs_tol", STRUMPACK_set_hss_abs_tol(*S,atol)); 180*34f43fa5SPieter Ghysels PetscFunctionReturn(0); 181*34f43fa5SPieter Ghysels } 182*34f43fa5SPieter Ghysels #undef __FUNCT__ 183*34f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSAbsTol" 184*34f43fa5SPieter Ghysels /*@ 185*34f43fa5SPieter Ghysels MatSTRUMPACKSetHSSAbsTol - Set STRUMPACK absolute tolerance for HSS compression 186*34f43fa5SPieter Ghysels Logically Collective on Mat 187*34f43fa5SPieter Ghysels 188*34f43fa5SPieter Ghysels Input Parameters: 189*34f43fa5SPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 190*34f43fa5SPieter Ghysels - atol - absolute compression tolerance 191*34f43fa5SPieter Ghysels 192*34f43fa5SPieter Ghysels Options Database: 193*34f43fa5SPieter Ghysels . -mat_strumpack_hss_abs_tol <1e-8> - Absolute compression tolerance (None) 194*34f43fa5SPieter Ghysels 195*34f43fa5SPieter Ghysels Level: beginner 196*34f43fa5SPieter Ghysels 197*34f43fa5SPieter Ghysels References: 198*34f43fa5SPieter Ghysels . STRUMPACK manual 199*34f43fa5SPieter Ghysels 200*34f43fa5SPieter Ghysels .seealso: MatGetFactor() 201*34f43fa5SPieter Ghysels @*/ 202*34f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSAbsTol(Mat F,PetscReal atol) 203*34f43fa5SPieter Ghysels { 204*34f43fa5SPieter Ghysels PetscErrorCode ierr; 205*34f43fa5SPieter Ghysels 206*34f43fa5SPieter Ghysels PetscFunctionBegin; 207*34f43fa5SPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 208*34f43fa5SPieter Ghysels PetscValidLogicalCollectiveReal(F,atol,2); 209*34f43fa5SPieter Ghysels ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSAbsTol_C",(Mat,PetscReal),(F,atol));CHKERRQ(ierr); 210*34f43fa5SPieter Ghysels PetscFunctionReturn(0); 211*34f43fa5SPieter Ghysels } 212*34f43fa5SPieter Ghysels 213*34f43fa5SPieter Ghysels 214*34f43fa5SPieter Ghysels #undef __FUNCT__ 215*34f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSMaxRank_STRUMPACK" 216*34f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSMaxRank_STRUMPACK(Mat F,PetscInt hssmaxrank) 217*34f43fa5SPieter Ghysels { 218*34f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 219*34f43fa5SPieter Ghysels 220*34f43fa5SPieter Ghysels PetscFunctionBegin; 221*34f43fa5SPieter Ghysels PetscStackCall("STRUMPACK_set_max_rank", STRUMPACK_set_max_rank(*S,hssmaxrank)); 222*34f43fa5SPieter Ghysels PetscFunctionReturn(0); 223*34f43fa5SPieter Ghysels } 224*34f43fa5SPieter Ghysels #undef __FUNCT__ 225*34f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSMaxRank" 226*34f43fa5SPieter Ghysels /*@ 227*34f43fa5SPieter Ghysels MatSTRUMPACKSetHSSMaxRank - Set STRUMPACK maximum HSS rank 228*34f43fa5SPieter Ghysels Logically Collective on Mat 229*34f43fa5SPieter Ghysels 230*34f43fa5SPieter Ghysels Input Parameters: 231*34f43fa5SPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 232*34f43fa5SPieter Ghysels - hssmaxrank - maximum rank used in low-rank approximation 233*34f43fa5SPieter Ghysels 234*34f43fa5SPieter Ghysels Options Database: 235*34f43fa5SPieter Ghysels . -mat_strumpack_max_rank - Maximum rank in HSS compression, when using pctype ilu (None) 236*34f43fa5SPieter Ghysels 237*34f43fa5SPieter Ghysels Level: beginner 238*34f43fa5SPieter Ghysels 239*34f43fa5SPieter Ghysels References: 240*34f43fa5SPieter Ghysels . STRUMPACK manual 241*34f43fa5SPieter Ghysels 242*34f43fa5SPieter Ghysels .seealso: MatGetFactor() 243*34f43fa5SPieter Ghysels @*/ 244*34f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSMaxRank(Mat F,PetscInt hssmaxrank) 245*34f43fa5SPieter Ghysels { 246*34f43fa5SPieter Ghysels PetscErrorCode ierr; 247*34f43fa5SPieter Ghysels 248*34f43fa5SPieter Ghysels PetscFunctionBegin; 249*34f43fa5SPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 250*34f43fa5SPieter Ghysels PetscValidLogicalCollectiveInt(F,hssmaxrank,2); 251*34f43fa5SPieter Ghysels ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSMaxRank_C",(Mat,PetscInt),(F,hssmaxrank));CHKERRQ(ierr); 252*34f43fa5SPieter Ghysels PetscFunctionReturn(0); 253*34f43fa5SPieter Ghysels } 254*34f43fa5SPieter Ghysels 255*34f43fa5SPieter Ghysels 256*34f43fa5SPieter Ghysels #undef __FUNCT__ 257291d0ed5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSMinFrontSize_STRUMPACK" 258291d0ed5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSMinFrontSize_STRUMPACK(Mat F,PetscInt hssminsize) 259575704cbSPieter Ghysels { 2600d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 261575704cbSPieter Ghysels 262575704cbSPieter Ghysels PetscFunctionBegin; 263291d0ed5SPieter Ghysels PetscStackCall("STRUMPACK_set_HSS_min_front_size", STRUMPACK_set_HSS_min_front_size(*S,hssminsize)); 264291d0ed5SPieter Ghysels PetscFunctionReturn(0); 265291d0ed5SPieter Ghysels } 266291d0ed5SPieter Ghysels #undef __FUNCT__ 267291d0ed5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSMinFrontSize" 268575704cbSPieter Ghysels /*@ 269291d0ed5SPieter Ghysels MatSTRUMPACKSetHSSMinFrontSize - Set STRUMPACK minimum dense matrix size for low-rank approximation 270575704cbSPieter Ghysels Logically Collective on Mat 271575704cbSPieter Ghysels 272575704cbSPieter Ghysels Input Parameters: 273575704cbSPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 274575704cbSPieter Ghysels - hssminsize - minimum dense matrix size for low-rank approximation 275575704cbSPieter Ghysels 276575704cbSPieter Ghysels Options Database: 277291d0ed5SPieter Ghysels . -mat_strumpack_hss_min_front_size <hssminsize> 278575704cbSPieter Ghysels 279575704cbSPieter Ghysels Level: beginner 280575704cbSPieter Ghysels 281575704cbSPieter Ghysels References: 282575704cbSPieter Ghysels . STRUMPACK manual 283575704cbSPieter Ghysels 284575704cbSPieter Ghysels .seealso: MatGetFactor() 285575704cbSPieter Ghysels @*/ 286291d0ed5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSMinFrontSize(Mat F,PetscInt hssminsize) 287575704cbSPieter Ghysels { 288575704cbSPieter Ghysels PetscErrorCode ierr; 289575704cbSPieter Ghysels 290575704cbSPieter Ghysels PetscFunctionBegin; 291575704cbSPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 292575704cbSPieter Ghysels PetscValidLogicalCollectiveInt(F,hssminsize,2); 293291d0ed5SPieter Ghysels ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSMinFrontSize_C",(Mat,PetscInt),(F,hssminsize));CHKERRQ(ierr); 294575704cbSPieter Ghysels PetscFunctionReturn(0); 295575704cbSPieter Ghysels } 296575704cbSPieter Ghysels 297291d0ed5SPieter Ghysels #undef __FUNCT__ 298*34f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSMinSepSize_STRUMPACK" 299*34f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSMinSepSize_STRUMPACK(Mat F,PetscInt hssminsize) 300*34f43fa5SPieter Ghysels { 301*34f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 302*34f43fa5SPieter Ghysels 303*34f43fa5SPieter Ghysels PetscFunctionBegin; 304*34f43fa5SPieter Ghysels PetscStackCall("STRUMPACK_set_HSS_min_sep_size", STRUMPACK_set_HSS_min_sep_size(*S,hssminsize)); 305*34f43fa5SPieter Ghysels PetscFunctionReturn(0); 306*34f43fa5SPieter Ghysels } 307*34f43fa5SPieter Ghysels #undef __FUNCT__ 308291d0ed5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSMinSepSize" 309291d0ed5SPieter Ghysels /*@ 310291d0ed5SPieter Ghysels MatSTRUMPACKSetHSSMinSepSize - Set STRUMPACK minimum separator size for low-rank approximation 311291d0ed5SPieter Ghysels Logically Collective on Mat 312291d0ed5SPieter Ghysels 313291d0ed5SPieter Ghysels Input Parameters: 314291d0ed5SPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 315291d0ed5SPieter Ghysels - hssminsize - minimum dense matrix size for low-rank approximation 316291d0ed5SPieter Ghysels 317291d0ed5SPieter Ghysels Options Database: 318291d0ed5SPieter Ghysels . -mat_strumpack_hss_min_sep_size <hssminsize> 319291d0ed5SPieter Ghysels 320291d0ed5SPieter Ghysels Level: beginner 321291d0ed5SPieter Ghysels 322291d0ed5SPieter Ghysels References: 323291d0ed5SPieter Ghysels . STRUMPACK manual 324291d0ed5SPieter Ghysels 325291d0ed5SPieter Ghysels .seealso: MatGetFactor() 326291d0ed5SPieter Ghysels @*/ 327291d0ed5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat F,PetscInt hssminsize) 328291d0ed5SPieter Ghysels { 329291d0ed5SPieter Ghysels PetscErrorCode ierr; 330291d0ed5SPieter Ghysels 331291d0ed5SPieter Ghysels PetscFunctionBegin; 332291d0ed5SPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 333291d0ed5SPieter Ghysels PetscValidLogicalCollectiveInt(F,hssminsize,2); 334291d0ed5SPieter Ghysels ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSMinSepSize_C",(Mat,PetscInt),(F,hssminsize));CHKERRQ(ierr); 335291d0ed5SPieter Ghysels PetscFunctionReturn(0); 336291d0ed5SPieter Ghysels } 337291d0ed5SPieter Ghysels 338291d0ed5SPieter Ghysels 339291d0ed5SPieter Ghysels #undef __FUNCT__ 340291d0ed5SPieter Ghysels #define __FUNCT__ "MatSolve_STRUMPACK" 341ad0c5e61SPieter Ghysels static PetscErrorCode MatSolve_STRUMPACK(Mat A,Vec b_mpi,Vec x) 3427d6ea485SPieter Ghysels { 3430d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)A->spptr; 3447d6ea485SPieter Ghysels STRUMPACK_RETURN_CODE sp_err; 3457d6ea485SPieter Ghysels PetscErrorCode ierr; 3467d6ea485SPieter Ghysels const PetscScalar *bptr; 3477d6ea485SPieter Ghysels PetscScalar *xptr; 3487d6ea485SPieter Ghysels 3497d6ea485SPieter Ghysels PetscFunctionBegin; 3507d6ea485SPieter Ghysels ierr = VecGetArray(x,&xptr);CHKERRQ(ierr); 3517d6ea485SPieter Ghysels ierr = VecGetArrayRead(b_mpi,&bptr);CHKERRQ(ierr); 3520d08a34dSPieter Ghysels 3530d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_solve",sp_err = STRUMPACK_solve(*S,(PetscScalar*)bptr,xptr,0)); 3540d08a34dSPieter Ghysels switch (sp_err) { 3550d08a34dSPieter Ghysels case STRUMPACK_SUCCESS: break; 3560d08a34dSPieter Ghysels case STRUMPACK_MATRIX_NOT_SET: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); break; } 3570d08a34dSPieter Ghysels case STRUMPACK_REORDERING_ERROR: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); break; } 3580d08a34dSPieter Ghysels default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: solve failed"); 3597d6ea485SPieter Ghysels } 3607d6ea485SPieter Ghysels ierr = VecRestoreArray(x,&xptr);CHKERRQ(ierr); 3617d6ea485SPieter Ghysels ierr = VecRestoreArrayRead(b_mpi,&bptr);CHKERRQ(ierr); 3627d6ea485SPieter Ghysels PetscFunctionReturn(0); 3637d6ea485SPieter Ghysels } 3647d6ea485SPieter Ghysels 365291d0ed5SPieter Ghysels #undef __FUNCT__ 366291d0ed5SPieter Ghysels #define __FUNCT__ "MatMatSolve_STRUMPACK" 367ad0c5e61SPieter Ghysels static PetscErrorCode MatMatSolve_STRUMPACK(Mat A,Mat B_mpi,Mat X) 3687d6ea485SPieter Ghysels { 3697d6ea485SPieter Ghysels PetscErrorCode ierr; 3707d6ea485SPieter Ghysels PetscBool flg; 3717d6ea485SPieter Ghysels 3727d6ea485SPieter Ghysels PetscFunctionBegin; 3737d6ea485SPieter Ghysels ierr = PetscObjectTypeCompareAny((PetscObject)B_mpi,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 3747d6ea485SPieter Ghysels if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 3757d6ea485SPieter Ghysels ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 3767d6ea485SPieter Ghysels if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); 3777d6ea485SPieter Ghysels SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_STRUMPACK() is not implemented yet"); 3787d6ea485SPieter Ghysels PetscFunctionReturn(0); 3797d6ea485SPieter Ghysels } 3807d6ea485SPieter Ghysels 381291d0ed5SPieter Ghysels #undef __FUNCT__ 382291d0ed5SPieter Ghysels #define __FUNCT__ "MatFactorInfo_STRUMPACK" 383ad0c5e61SPieter Ghysels static PetscErrorCode MatFactorInfo_STRUMPACK(Mat A,PetscViewer viewer) 384ad0c5e61SPieter Ghysels { 385ad0c5e61SPieter Ghysels PetscErrorCode ierr; 386ad0c5e61SPieter Ghysels 387ad0c5e61SPieter Ghysels PetscFunctionBegin; 388ad0c5e61SPieter Ghysels /* check if matrix is strumpack type */ 389ad0c5e61SPieter Ghysels if (A->ops->solve != MatSolve_STRUMPACK) PetscFunctionReturn(0); 390ad0c5e61SPieter Ghysels ierr = PetscViewerASCIIPrintf(viewer,"STRUMPACK sparse solver!\n");CHKERRQ(ierr); 391ad0c5e61SPieter Ghysels PetscFunctionReturn(0); 392ad0c5e61SPieter Ghysels } 393ad0c5e61SPieter Ghysels 394291d0ed5SPieter Ghysels #undef __FUNCT__ 395291d0ed5SPieter Ghysels #define __FUNCT__ "MatView_STRUMPACK" 396ad0c5e61SPieter Ghysels static PetscErrorCode MatView_STRUMPACK(Mat A,PetscViewer viewer) 397ad0c5e61SPieter Ghysels { 398ad0c5e61SPieter Ghysels PetscErrorCode ierr; 399ad0c5e61SPieter Ghysels PetscBool iascii; 400ad0c5e61SPieter Ghysels PetscViewerFormat format; 401ad0c5e61SPieter Ghysels 402ad0c5e61SPieter Ghysels PetscFunctionBegin; 403ad0c5e61SPieter Ghysels ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 404ad0c5e61SPieter Ghysels if (iascii) { 405ad0c5e61SPieter Ghysels ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 406ad0c5e61SPieter Ghysels if (format == PETSC_VIEWER_ASCII_INFO) { 407ad0c5e61SPieter Ghysels ierr = MatFactorInfo_STRUMPACK(A,viewer);CHKERRQ(ierr); 408ad0c5e61SPieter Ghysels } 409ad0c5e61SPieter Ghysels } 410ad0c5e61SPieter Ghysels PetscFunctionReturn(0); 411ad0c5e61SPieter Ghysels } 4127d6ea485SPieter Ghysels 413291d0ed5SPieter Ghysels #undef __FUNCT__ 414291d0ed5SPieter Ghysels #define __FUNCT__ "MatLUFactorNumeric_STRUMPACK" 415ad0c5e61SPieter Ghysels static PetscErrorCode MatLUFactorNumeric_STRUMPACK(Mat F,Mat A,const MatFactorInfo *info) 4167d6ea485SPieter Ghysels { 4170d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 4187d6ea485SPieter Ghysels STRUMPACK_RETURN_CODE sp_err; 4197d6ea485SPieter Ghysels Mat_SeqAIJ *A_d,*A_o; 4207d6ea485SPieter Ghysels Mat_MPIAIJ *mat; 4217d6ea485SPieter Ghysels PetscErrorCode ierr; 4220d08a34dSPieter Ghysels PetscInt M=A->rmap->N,m=A->rmap->n; 4237d6ea485SPieter Ghysels PetscBool flg; 4247d6ea485SPieter Ghysels 4257d6ea485SPieter Ghysels PetscFunctionBegin; 4267d6ea485SPieter Ghysels ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); 4277d6ea485SPieter Ghysels if (flg) { /* A is MATMPIAIJ */ 4287d6ea485SPieter Ghysels mat = (Mat_MPIAIJ*)A->data; 4297d6ea485SPieter Ghysels A_d = (Mat_SeqAIJ*)(mat->A)->data; 4307d6ea485SPieter Ghysels A_o = (Mat_SeqAIJ*)(mat->B)->data; 4310d08a34dSPieter 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)); 4327d6ea485SPieter Ghysels } else { /* A is MATSEQAIJ */ 4337d6ea485SPieter Ghysels A_d = (Mat_SeqAIJ*)A->data; 4340d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_set_csr_matrix",STRUMPACK_set_csr_matrix(*S,&M,A_d->i,A_d->j,A_d->a,0)); 4357d6ea485SPieter Ghysels } 4367d6ea485SPieter Ghysels 4377d6ea485SPieter Ghysels /* Reorder and Factor the matrix. */ 4387d6ea485SPieter Ghysels /* TODO figure out how to avoid reorder if the matrix values changed, but the pattern remains the same. */ 4390d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_reorder",sp_err = STRUMPACK_reorder(*S)); 4400d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_factor",sp_err = STRUMPACK_factor(*S)); 4410d08a34dSPieter Ghysels switch (sp_err) { 4420d08a34dSPieter Ghysels case STRUMPACK_SUCCESS: break; 4430d08a34dSPieter Ghysels case STRUMPACK_MATRIX_NOT_SET: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); break; } 4440d08a34dSPieter Ghysels case STRUMPACK_REORDERING_ERROR: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); break; } 4450d08a34dSPieter Ghysels default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: factorization failed"); 4467d6ea485SPieter Ghysels } 4477d6ea485SPieter Ghysels PetscFunctionReturn(0); 4487d6ea485SPieter Ghysels } 4497d6ea485SPieter Ghysels 450291d0ed5SPieter Ghysels #undef __FUNCT__ 451291d0ed5SPieter Ghysels #define __FUNCT__ "MatLUFactorSymbolic_STRUMPACK" 452ad0c5e61SPieter Ghysels static PetscErrorCode MatLUFactorSymbolic_STRUMPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 4537d6ea485SPieter Ghysels { 4547d6ea485SPieter Ghysels PetscFunctionBegin; 4557d6ea485SPieter Ghysels F->ops->lufactornumeric = MatLUFactorNumeric_STRUMPACK; 4567d6ea485SPieter Ghysels F->ops->solve = MatSolve_STRUMPACK; 4577d6ea485SPieter Ghysels F->ops->matsolve = MatMatSolve_STRUMPACK; 4587d6ea485SPieter Ghysels PetscFunctionReturn(0); 4597d6ea485SPieter Ghysels } 4607d6ea485SPieter Ghysels 461291d0ed5SPieter Ghysels #undef __FUNCT__ 462291d0ed5SPieter Ghysels #define __FUNCT__ "MatFactorGetSolverPackage_aij_strumpack" 463ad0c5e61SPieter Ghysels static PetscErrorCode MatFactorGetSolverPackage_aij_strumpack(Mat A,const MatSolverPackage *type) 4647d6ea485SPieter Ghysels { 4657d6ea485SPieter Ghysels PetscFunctionBegin; 4667d6ea485SPieter Ghysels *type = MATSOLVERSTRUMPACK; 4677d6ea485SPieter Ghysels PetscFunctionReturn(0); 4687d6ea485SPieter Ghysels } 4697d6ea485SPieter Ghysels 470575704cbSPieter Ghysels /*MC 471575704cbSPieter Ghysels MATSOLVERSSTRUMPACK = "strumpack" - A solver package providing a direct sparse solver (PCLU) 472575704cbSPieter Ghysels and a preconditioner (PCILU) using low-rank compression via the external package STRUMPACK. 473575704cbSPieter Ghysels 474575704cbSPieter Ghysels Consult the STRUMPACK-sparse manual for more info. 475575704cbSPieter Ghysels 476575704cbSPieter Ghysels Use 477575704cbSPieter Ghysels ./configure --download-strumpack 478575704cbSPieter Ghysels to have PETSc installed with STRUMPACK 479575704cbSPieter Ghysels 480575704cbSPieter Ghysels Use 481575704cbSPieter Ghysels -pc_type lu -pc_factor_mat_solver_package strumpack 482575704cbSPieter Ghysels to use this as an exact (direct) solver, use 483575704cbSPieter Ghysels -pc_type ilu -pc_factor_mat_solver_package strumpack 484575704cbSPieter Ghysels to enable low-rank compression (i.e, use as a preconditioner). 485575704cbSPieter Ghysels 486575704cbSPieter Ghysels Works with AIJ matrices 487575704cbSPieter Ghysels 488575704cbSPieter Ghysels Options Database Keys: 489*34f43fa5SPieter Ghysels + -mat_strumpack_verbose 490*34f43fa5SPieter Ghysels . -mat_strumpack_hss_rel_tol <1e-2> - Relative compression tolerance (None) 491*34f43fa5SPieter Ghysels . -mat_strumpack_hss_abs_tol <1e-8> - Absolute compression tolerance (None) 492575704cbSPieter Ghysels . -mat_strumpack_colperm <TRUE> - Permute matrix to make diagonal nonzeros (None) 493*34f43fa5SPieter Ghysels . -mat_strumpack_hss_min_front_size <1000> - Minimum size of dense block for HSS compression (None) 494*34f43fa5SPieter Ghysels . -mat_strumpack_hss_min_sep_size <256> - Minimum size of separator for HSS compression (None) 495*34f43fa5SPieter Ghysels . -mat_strumpack_max_rank - Maximum rank in HSS compression, when using pctype ilu (None) 496*34f43fa5SPieter Ghysels . -mat_strumpack_reordering <METIS> - Sparsity reducing matrix reordering (choose one of) NATURAL METIS PARMETIS SCOTCH PTSCOTCH RCM (None) 497*34f43fa5SPieter Ghysels . -mat_strumpack_iterative_solver <DIRECT> - Select iterative solver from STRUMPACK (choose one of) AUTO DIRECT REFINE PREC_GMRES GMRES PREC_BICGSTAB BICGSTAB (None) 498575704cbSPieter Ghysels 499575704cbSPieter Ghysels Level: beginner 500575704cbSPieter Ghysels 501575704cbSPieter Ghysels .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, PCFactorSetMatSolverPackage(), MatSolverPackage 502575704cbSPieter Ghysels M*/ 503291d0ed5SPieter Ghysels #undef __FUNCT__ 504291d0ed5SPieter Ghysels #define __FUNCT__ "MatGetFactor_aij_strumpack" 505ad0c5e61SPieter Ghysels static PetscErrorCode MatGetFactor_aij_strumpack(Mat A,MatFactorType ftype,Mat *F) 5067d6ea485SPieter Ghysels { 5077d6ea485SPieter Ghysels Mat B; 5087d6ea485SPieter Ghysels PetscErrorCode ierr; 5097d6ea485SPieter Ghysels PetscInt M=A->rmap->N,N=A->cmap->N; 510575704cbSPieter Ghysels PetscBool verb,flg,set; 511*34f43fa5SPieter Ghysels PetscReal ctol; 512*34f43fa5SPieter Ghysels PetscInt hssminsize,max_rank; 513*34f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S; 514*34f43fa5SPieter Ghysels STRUMPACK_INTERFACE iface; 515*34f43fa5SPieter Ghysels STRUMPACK_REORDERING_STRATEGY ndcurrent,ndvalue; 516*34f43fa5SPieter Ghysels STRUMPACK_KRYLOV_SOLVER itcurrent,itsolver; 517*34f43fa5SPieter Ghysels const STRUMPACK_PRECISION table[2][2][2] = 518*34f43fa5SPieter Ghysels {{{STRUMPACK_FLOATCOMPLEX_64, STRUMPACK_DOUBLECOMPLEX_64}, 5197d6ea485SPieter Ghysels {STRUMPACK_FLOAT_64, STRUMPACK_DOUBLE_64}}, 5207d6ea485SPieter Ghysels {{STRUMPACK_FLOATCOMPLEX, STRUMPACK_DOUBLECOMPLEX}, 5217d6ea485SPieter Ghysels {STRUMPACK_FLOAT, STRUMPACK_DOUBLE}}}; 522*34f43fa5SPieter Ghysels const STRUMPACK_PRECISION prec = 523*34f43fa5SPieter Ghysels table[(sizeof(PetscInt)==8)?0:1] 524*34f43fa5SPieter Ghysels [(PETSC_SCALAR==PETSC_COMPLEX)?0:1] 525*34f43fa5SPieter Ghysels [(PETSC_REAL==PETSC_FLOAT)?0:1]; 5267d6ea485SPieter Ghysels 5277d6ea485SPieter Ghysels PetscFunctionBegin; 5287d6ea485SPieter Ghysels /* Create the factorization matrix */ 5297d6ea485SPieter Ghysels ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 5307d6ea485SPieter Ghysels ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,M,N);CHKERRQ(ierr); 5317d6ea485SPieter Ghysels ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 5327d6ea485SPieter Ghysels ierr = MatSeqAIJSetPreallocation(B,0,NULL); 5337d6ea485SPieter Ghysels ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr); 534575704cbSPieter Ghysels if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) { 5357d6ea485SPieter Ghysels B->ops->lufactorsymbolic = MatLUFactorSymbolic_STRUMPACK; 536575704cbSPieter Ghysels B->ops->ilufactorsymbolic = MatLUFactorSymbolic_STRUMPACK; 537575704cbSPieter Ghysels } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 5387d6ea485SPieter Ghysels B->ops->view = MatView_STRUMPACK; 5397d6ea485SPieter Ghysels B->ops->destroy = MatDestroy_STRUMPACK; 5407d6ea485SPieter Ghysels B->ops->getdiagonal = MatGetDiagonal_STRUMPACK; 5417d6ea485SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_aij_strumpack);CHKERRQ(ierr); 542*34f43fa5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetReordering_C",MatSTRUMPACKSetReordering_STRUMPACK);CHKERRQ(ierr); 543*34f43fa5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetColPerm_C",MatSTRUMPACKSetColPerm_STRUMPACK);CHKERRQ(ierr); 544*34f43fa5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSRelTol_C",MatSTRUMPACKSetHSSRelTol_STRUMPACK);CHKERRQ(ierr); 545*34f43fa5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSAbsTol_C",MatSTRUMPACKSetHSSAbsTol_STRUMPACK);CHKERRQ(ierr); 546291d0ed5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSMinFrontSize_C",MatSTRUMPACKSetHSSMinFrontSize_STRUMPACK);CHKERRQ(ierr); 547291d0ed5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSMinSepSize_C",MatSTRUMPACKSetHSSMinSepSize_STRUMPACK);CHKERRQ(ierr); 548*34f43fa5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSMaxRank_C",MatSTRUMPACKSetHSSMaxRank_STRUMPACK);CHKERRQ(ierr); 549575704cbSPieter Ghysels B->factortype = ftype; 5500d08a34dSPieter Ghysels ierr = PetscNewLog(B,&S);CHKERRQ(ierr); 5510d08a34dSPieter Ghysels B->spptr = S; 5520d08a34dSPieter Ghysels 5530d08a34dSPieter Ghysels ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg); 5540d08a34dSPieter Ghysels iface = flg ? STRUMPACK_MT : STRUMPACK_MPI_DIST; 5557d6ea485SPieter Ghysels 5567d6ea485SPieter Ghysels ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"STRUMPACK Options","Mat");CHKERRQ(ierr); 5577d6ea485SPieter Ghysels 558575704cbSPieter Ghysels verb = PetscLogPrintInfo ? PETSC_TRUE : PETSC_FALSE; 5597d6ea485SPieter Ghysels ierr = PetscOptionsBool("-mat_strumpack_verbose","Print STRUMPACK information","None",verb,&verb,NULL);CHKERRQ(ierr); 5607d6ea485SPieter Ghysels 561*34f43fa5SPieter Ghysels PetscStackCall("STRUMPACK_init",STRUMPACK_init(S,PetscObjectComm((PetscObject)A),prec,iface,0,NULL,verb)); 56255c022e5SPieter Ghysels 563*34f43fa5SPieter Ghysels PetscStackCall("STRUMPACK_hss_rel_tol",ctol = (PetscReal)STRUMPACK_hss_rel_tol(*S)); 564*34f43fa5SPieter Ghysels ierr = PetscOptionsReal("-mat_strumpack_hss_rel_tol","Relative compression tolerance","None",ctol,&ctol,&set);CHKERRQ(ierr); 565*34f43fa5SPieter Ghysels if (set) PetscStackCall("STRUMPACK_set_hss_rel_tol",STRUMPACK_set_hss_rel_tol(*S,(double)ctol)); 5667d6ea485SPieter Ghysels 567*34f43fa5SPieter Ghysels PetscStackCall("STRUMPACK_hss_abs_tol",ctol = (PetscReal)STRUMPACK_hss_abs_tol(*S)); 568*34f43fa5SPieter Ghysels ierr = PetscOptionsReal("-mat_strumpack_hss_abs_tol","Absolute compression tolerance","None",ctol,&ctol,&set);CHKERRQ(ierr); 569*34f43fa5SPieter Ghysels if (set) PetscStackCall("STRUMPACK_set_hss_abs_tol",STRUMPACK_set_hss_abs_tol(*S,(double)ctol)); 570575704cbSPieter Ghysels 571291d0ed5SPieter Ghysels PetscStackCall("STRUMPACK_mc64job",flg = (STRUMPACK_mc64job(*S) == 0) ? PETSC_FALSE : PETSC_TRUE); 572575704cbSPieter Ghysels ierr = PetscOptionsBool("-mat_strumpack_colperm","Find a col perm to get nonzero diagonal","None",flg,&flg,&set);CHKERRQ(ierr); 5730d08a34dSPieter Ghysels if (set) PetscStackCall("STRUMPACK_set_mc64job",STRUMPACK_set_mc64job(*S,flg ? 5 : 0)); 574575704cbSPieter Ghysels 575291d0ed5SPieter Ghysels PetscStackCall("STRUMPACK_HSS_min_front_size",hssminsize = (PetscInt)STRUMPACK_HSS_min_front_size(*S)); 576291d0ed5SPieter Ghysels ierr = PetscOptionsInt("-mat_strumpack_hss_min_front_size","Minimum size of dense block for HSS compression","None",hssminsize,&hssminsize,&set);CHKERRQ(ierr); 577291d0ed5SPieter Ghysels if (set) PetscStackCall("STRUMPACK_set_HSS_min_front_size",STRUMPACK_set_HSS_min_front_size(*S,(int)hssminsize)); 578291d0ed5SPieter Ghysels 579291d0ed5SPieter Ghysels PetscStackCall("STRUMPACK_HSS_min_sep_size",hssminsize = (PetscInt)STRUMPACK_HSS_min_sep_size(*S)); 580291d0ed5SPieter Ghysels ierr = PetscOptionsInt("-mat_strumpack_hss_min_sep_size","Minimum size of separator for HSS compression","None",hssminsize,&hssminsize,&set);CHKERRQ(ierr); 581291d0ed5SPieter Ghysels if (set) PetscStackCall("STRUMPACK_set_HSS_min_sep_size",STRUMPACK_set_HSS_min_sep_size(*S,(int)hssminsize)); 582575704cbSPieter Ghysels 583*34f43fa5SPieter Ghysels PetscStackCall("STRUMPACK_max_rank",max_rank = (PetscInt)STRUMPACK_max_rank(*S)); 584*34f43fa5SPieter Ghysels ierr = PetscOptionsInt("-mat_strumpack_max_rank","Maximum rank in HSS approximation","None",max_rank,&max_rank,&set);CHKERRQ(ierr); 585*34f43fa5SPieter Ghysels if (set) PetscStackCall("STRUMPACK_set_max_rank",STRUMPACK_set_max_rank(*S,(int)max_rank)); 586*34f43fa5SPieter Ghysels 587*34f43fa5SPieter Ghysels const char *const STRUMPACKNDTypes[] = {"NATURAL","METIS","PARMETIS","SCOTCH","PTSCOTCH","RCM","STRUMPACKNDTypes","",0}; 588*34f43fa5SPieter Ghysels PetscStackCall("STRUMPACK_reordering_method",ndcurrent = STRUMPACK_reordering_method(*S)); 589*34f43fa5SPieter Ghysels PetscOptionsEnum("-mat_strumpack_reordering","Sparsity reducing matrix reordering","None",STRUMPACKNDTypes,(PetscEnum)ndcurrent,(PetscEnum*)&ndvalue,&set);CHKERRQ(ierr); 590*34f43fa5SPieter Ghysels if (set) PetscStackCall("STRUMPACK_set_reordering_method",STRUMPACK_set_reordering_method(*S,ndvalue)); 591575704cbSPieter Ghysels 59288382b05SPieter Ghysels if (ftype == MAT_FACTOR_ILU) { 59388382b05SPieter Ghysels /* When enabling HSS compression, the STRUMPACK solver becomes an incomplete */ 59488382b05SPieter Ghysels /* (or approximate) LU factorization. */ 595291d0ed5SPieter Ghysels PetscStackCall("STRUMPACK_enable_HSS",STRUMPACK_enable_HSS(*S)); 59688382b05SPieter Ghysels } 597575704cbSPieter Ghysels 598*34f43fa5SPieter Ghysels /* Disable the outer iterative solver from STRUMPACK. */ 599*34f43fa5SPieter Ghysels /* When STRUMPACK is used as a direct solver, it will by default do iterative refinement. */ 600*34f43fa5SPieter Ghysels /* When STRUMPACK is used as an approximate factorization preconditioner (by enabling */ 601*34f43fa5SPieter Ghysels /* low-rank compression), it will use it's own preconditioned GMRES. Here we can disable */ 602*34f43fa5SPieter Ghysels /* the outer iterative solver, as PETSc uses STRUMPACK from within a KSP. */ 603*34f43fa5SPieter Ghysels PetscStackCall("STRUMPACK_set_Krylov_solver", STRUMPACK_set_Krylov_solver(*S, STRUMPACK_DIRECT)); 604*34f43fa5SPieter Ghysels 605*34f43fa5SPieter Ghysels const char *const SolverTypes[] = {"AUTO","NONE","REFINE","PREC_GMRES","GMRES","PREC_BICGSTAB","BICGSTAB","SolverTypes","",0}; 606*34f43fa5SPieter Ghysels PetscStackCall("STRUMPACK_Krylov_solver",itcurrent = STRUMPACK_Krylov_solver(*S)); 607*34f43fa5SPieter Ghysels PetscOptionsEnum("-mat_strumpack_iterative_solver","Select iterative solver from STRUMPACK","None",SolverTypes,(PetscEnum)itcurrent,(PetscEnum*)&itsolver,&set);CHKERRQ(ierr); 608*34f43fa5SPieter Ghysels if (set) PetscStackCall("STRUMPACK_set_Krylov_solver",STRUMPACK_set_Krylov_solver(*S,itsolver)); 609*34f43fa5SPieter Ghysels 610*34f43fa5SPieter Ghysels PetscOptionsEnd(); 61155c022e5SPieter Ghysels 6127d6ea485SPieter Ghysels *F = B; 6137d6ea485SPieter Ghysels PetscFunctionReturn(0); 6147d6ea485SPieter Ghysels } 6157d6ea485SPieter Ghysels 616291d0ed5SPieter Ghysels #undef __FUNCT__ 617291d0ed5SPieter Ghysels #define __FUNCT__ "MatSolverPackageRegister_STRUMPACK" 6187d6ea485SPieter Ghysels PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_STRUMPACK(void) 6197d6ea485SPieter Ghysels { 6207d6ea485SPieter Ghysels PetscErrorCode ierr; 6217d6ea485SPieter Ghysels 6227d6ea485SPieter Ghysels PetscFunctionBegin; 6237d6ea485SPieter Ghysels ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr); 6247d6ea485SPieter Ghysels ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr); 625575704cbSPieter Ghysels ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_ILU,MatGetFactor_aij_strumpack);CHKERRQ(ierr); 626575704cbSPieter Ghysels ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_ILU,MatGetFactor_aij_strumpack);CHKERRQ(ierr); 6277d6ea485SPieter Ghysels PetscFunctionReturn(0); 6287d6ea485SPieter Ghysels } 629