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 5*291d0ed5SPieter Ghysels #undef __FUNCT__ 6*291d0ed5SPieter 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 14*291d0ed5SPieter Ghysels #undef __FUNCT__ 15*291d0ed5SPieter 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); 35575704cbSPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSRelCompTol_C",NULL);CHKERRQ(ierr); 36*291d0ed5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSMinFrontSize_C",NULL);CHKERRQ(ierr); 37*291d0ed5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSMinSepSize_C",NULL);CHKERRQ(ierr); 38575704cbSPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetColPerm_C",NULL);CHKERRQ(ierr); 39575704cbSPieter Ghysels 40575704cbSPieter Ghysels PetscFunctionReturn(0); 41575704cbSPieter Ghysels } 42575704cbSPieter Ghysels 43*291d0ed5SPieter Ghysels #undef __FUNCT__ 44*291d0ed5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSRelCompTol_STRUMPACK" 45575704cbSPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSRelCompTol_STRUMPACK(Mat F,PetscReal rctol) 46575704cbSPieter Ghysels { 470d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 48575704cbSPieter Ghysels 49575704cbSPieter Ghysels PetscFunctionBegin; 50*291d0ed5SPieter Ghysels PetscStackCall("STRUMPACK_set_hss_rel_tol", STRUMPACK_set_hss_rel_tol(*S,rctol)); 51575704cbSPieter Ghysels PetscFunctionReturn(0); 52575704cbSPieter Ghysels } 53575704cbSPieter Ghysels 54*291d0ed5SPieter Ghysels #undef __FUNCT__ 55*291d0ed5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSRelCompTol" 56575704cbSPieter Ghysels /*@ 57575704cbSPieter Ghysels MatSTRUMPACKSetHSSRelCompTol - Set STRUMPACK relative tolerance for HSS compression 58575704cbSPieter Ghysels Logically Collective on Mat 59575704cbSPieter Ghysels 60575704cbSPieter Ghysels Input Parameters: 61575704cbSPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 62575704cbSPieter Ghysels - rctol - relative compression tolerance 63575704cbSPieter Ghysels 64575704cbSPieter Ghysels Options Database: 65575704cbSPieter Ghysels . -mat_strumpack_rctol <rctol> 66575704cbSPieter Ghysels 67575704cbSPieter Ghysels Level: beginner 68575704cbSPieter Ghysels 69575704cbSPieter Ghysels References: 70575704cbSPieter Ghysels . STRUMPACK manual 71575704cbSPieter Ghysels 72575704cbSPieter Ghysels .seealso: MatGetFactor() 73575704cbSPieter Ghysels @*/ 74575704cbSPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat F,PetscReal rctol) 75575704cbSPieter Ghysels { 76575704cbSPieter Ghysels PetscErrorCode ierr; 77575704cbSPieter Ghysels 78575704cbSPieter Ghysels PetscFunctionBegin; 79575704cbSPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 80575704cbSPieter Ghysels PetscValidLogicalCollectiveReal(F,rctol,2); 81575704cbSPieter Ghysels ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSRelCompTol_C",(Mat,PetscReal),(F,rctol));CHKERRQ(ierr); 82575704cbSPieter Ghysels PetscFunctionReturn(0); 83575704cbSPieter Ghysels } 84575704cbSPieter Ghysels 85*291d0ed5SPieter Ghysels #undef __FUNCT__ 86*291d0ed5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSMinFrontSize_STRUMPACK" 87*291d0ed5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSMinFrontSize_STRUMPACK(Mat F,PetscInt hssminsize) 88575704cbSPieter Ghysels { 890d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 90575704cbSPieter Ghysels 91575704cbSPieter Ghysels PetscFunctionBegin; 92*291d0ed5SPieter Ghysels PetscStackCall("STRUMPACK_set_HSS_min_front_size", STRUMPACK_set_HSS_min_front_size(*S,hssminsize)); 93*291d0ed5SPieter Ghysels PetscFunctionReturn(0); 94*291d0ed5SPieter Ghysels } 95*291d0ed5SPieter Ghysels #undef __FUNCT__ 96*291d0ed5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSMinSepSize_STRUMPACK" 97*291d0ed5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSMinSepSize_STRUMPACK(Mat F,PetscInt hssminsize) 98*291d0ed5SPieter Ghysels { 99*291d0ed5SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 100*291d0ed5SPieter Ghysels 101*291d0ed5SPieter Ghysels PetscFunctionBegin; 102*291d0ed5SPieter Ghysels PetscStackCall("STRUMPACK_set_HSS_min_sep_size", STRUMPACK_set_HSS_min_sep_size(*S,hssminsize)); 103575704cbSPieter Ghysels PetscFunctionReturn(0); 104575704cbSPieter Ghysels } 105575704cbSPieter Ghysels 106*291d0ed5SPieter Ghysels #undef __FUNCT__ 107*291d0ed5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSMinFrontSize" 108575704cbSPieter Ghysels /*@ 109*291d0ed5SPieter Ghysels MatSTRUMPACKSetHSSMinFrontSize - Set STRUMPACK minimum dense matrix size for low-rank approximation 110575704cbSPieter Ghysels Logically Collective on Mat 111575704cbSPieter Ghysels 112575704cbSPieter Ghysels Input Parameters: 113575704cbSPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 114575704cbSPieter Ghysels - hssminsize - minimum dense matrix size for low-rank approximation 115575704cbSPieter Ghysels 116575704cbSPieter Ghysels Options Database: 117*291d0ed5SPieter Ghysels . -mat_strumpack_hss_min_front_size <hssminsize> 118575704cbSPieter Ghysels 119575704cbSPieter Ghysels Level: beginner 120575704cbSPieter Ghysels 121575704cbSPieter Ghysels References: 122575704cbSPieter Ghysels . STRUMPACK manual 123575704cbSPieter Ghysels 124575704cbSPieter Ghysels .seealso: MatGetFactor() 125575704cbSPieter Ghysels @*/ 126*291d0ed5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSMinFrontSize(Mat F,PetscInt hssminsize) 127575704cbSPieter Ghysels { 128575704cbSPieter Ghysels PetscErrorCode ierr; 129575704cbSPieter Ghysels 130575704cbSPieter Ghysels PetscFunctionBegin; 131575704cbSPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 132575704cbSPieter Ghysels PetscValidLogicalCollectiveInt(F,hssminsize,2); 133*291d0ed5SPieter Ghysels ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSMinFrontSize_C",(Mat,PetscInt),(F,hssminsize));CHKERRQ(ierr); 134575704cbSPieter Ghysels PetscFunctionReturn(0); 135575704cbSPieter Ghysels } 136575704cbSPieter Ghysels 137*291d0ed5SPieter Ghysels #undef __FUNCT__ 138*291d0ed5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSMinSepSize" 139*291d0ed5SPieter Ghysels /*@ 140*291d0ed5SPieter Ghysels MatSTRUMPACKSetHSSMinSepSize - Set STRUMPACK minimum separator size for low-rank approximation 141*291d0ed5SPieter Ghysels Logically Collective on Mat 142*291d0ed5SPieter Ghysels 143*291d0ed5SPieter Ghysels Input Parameters: 144*291d0ed5SPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 145*291d0ed5SPieter Ghysels - hssminsize - minimum dense matrix size for low-rank approximation 146*291d0ed5SPieter Ghysels 147*291d0ed5SPieter Ghysels Options Database: 148*291d0ed5SPieter Ghysels . -mat_strumpack_hss_min_sep_size <hssminsize> 149*291d0ed5SPieter Ghysels 150*291d0ed5SPieter Ghysels Level: beginner 151*291d0ed5SPieter Ghysels 152*291d0ed5SPieter Ghysels References: 153*291d0ed5SPieter Ghysels . STRUMPACK manual 154*291d0ed5SPieter Ghysels 155*291d0ed5SPieter Ghysels .seealso: MatGetFactor() 156*291d0ed5SPieter Ghysels @*/ 157*291d0ed5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat F,PetscInt hssminsize) 158*291d0ed5SPieter Ghysels { 159*291d0ed5SPieter Ghysels PetscErrorCode ierr; 160*291d0ed5SPieter Ghysels 161*291d0ed5SPieter Ghysels PetscFunctionBegin; 162*291d0ed5SPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 163*291d0ed5SPieter Ghysels PetscValidLogicalCollectiveInt(F,hssminsize,2); 164*291d0ed5SPieter Ghysels ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSMinSepSize_C",(Mat,PetscInt),(F,hssminsize));CHKERRQ(ierr); 165*291d0ed5SPieter Ghysels PetscFunctionReturn(0); 166*291d0ed5SPieter Ghysels } 167*291d0ed5SPieter Ghysels 168*291d0ed5SPieter Ghysels 169*291d0ed5SPieter Ghysels #undef __FUNCT__ 170*291d0ed5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetColPerm_STRUMPACK" 171575704cbSPieter Ghysels static PetscErrorCode MatSTRUMPACKSetColPerm_STRUMPACK(Mat F,PetscBool cperm) 172575704cbSPieter Ghysels { 1730d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 174575704cbSPieter Ghysels 175575704cbSPieter Ghysels PetscFunctionBegin; 1760d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_set_mc64job",STRUMPACK_set_mc64job(*S,cperm ? 5 : 0)); 177575704cbSPieter Ghysels PetscFunctionReturn(0); 178575704cbSPieter Ghysels } 179575704cbSPieter Ghysels 180*291d0ed5SPieter Ghysels #undef __FUNCT__ 181*291d0ed5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetColPerm" 182575704cbSPieter Ghysels /*@ 183575704cbSPieter Ghysels MatSTRUMPACKSetColPerm - Set whether STRUMPACK should try to permute the columns of the matrix in order to get a nonzero diagonal 184575704cbSPieter Ghysels Logically Collective on Mat 185575704cbSPieter Ghysels 186575704cbSPieter Ghysels Input Parameters: 187575704cbSPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 188575704cbSPieter Ghysels - cperm - whether or not to permute (internally) the columns of the matrix 189575704cbSPieter Ghysels 190575704cbSPieter Ghysels Options Database: 191575704cbSPieter Ghysels . -mat_strumpack_colperm <cperm> 192575704cbSPieter Ghysels 193575704cbSPieter Ghysels Level: beginner 194575704cbSPieter Ghysels 195575704cbSPieter Ghysels References: 196575704cbSPieter Ghysels . STRUMPACK manual 197575704cbSPieter Ghysels 198575704cbSPieter Ghysels .seealso: MatGetFactor() 199575704cbSPieter Ghysels @*/ 200575704cbSPieter Ghysels PetscErrorCode MatSTRUMPACKSetColPerm(Mat F,PetscBool cperm) 201575704cbSPieter Ghysels { 202575704cbSPieter Ghysels PetscErrorCode ierr; 203575704cbSPieter Ghysels 204575704cbSPieter Ghysels PetscFunctionBegin; 205575704cbSPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 206575704cbSPieter Ghysels PetscValidLogicalCollectiveBool(F,cperm,2); 207575704cbSPieter Ghysels ierr = PetscTryMethod(F,"MatSTRUMPACKSetColPerm_C",(Mat,PetscBool),(F,cperm));CHKERRQ(ierr); 2087d6ea485SPieter Ghysels PetscFunctionReturn(0); 2097d6ea485SPieter Ghysels } 2107d6ea485SPieter Ghysels 211*291d0ed5SPieter Ghysels #undef __FUNCT__ 212*291d0ed5SPieter Ghysels #define __FUNCT__ "MatSolve_STRUMPACK" 213ad0c5e61SPieter Ghysels static PetscErrorCode MatSolve_STRUMPACK(Mat A,Vec b_mpi,Vec x) 2147d6ea485SPieter Ghysels { 2150d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)A->spptr; 2167d6ea485SPieter Ghysels STRUMPACK_RETURN_CODE sp_err; 2177d6ea485SPieter Ghysels PetscErrorCode ierr; 2187d6ea485SPieter Ghysels const PetscScalar *bptr; 2197d6ea485SPieter Ghysels PetscScalar *xptr; 2207d6ea485SPieter Ghysels 2217d6ea485SPieter Ghysels PetscFunctionBegin; 2227d6ea485SPieter Ghysels ierr = VecGetArray(x,&xptr);CHKERRQ(ierr); 2237d6ea485SPieter Ghysels ierr = VecGetArrayRead(b_mpi,&bptr);CHKERRQ(ierr); 2240d08a34dSPieter Ghysels 2250d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_solve",sp_err = STRUMPACK_solve(*S,(PetscScalar*)bptr,xptr,0)); 2260d08a34dSPieter Ghysels switch (sp_err) { 2270d08a34dSPieter Ghysels case STRUMPACK_SUCCESS: break; 2280d08a34dSPieter Ghysels case STRUMPACK_MATRIX_NOT_SET: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); break; } 2290d08a34dSPieter Ghysels case STRUMPACK_REORDERING_ERROR: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); break; } 2300d08a34dSPieter Ghysels default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: solve failed"); 2317d6ea485SPieter Ghysels } 2327d6ea485SPieter Ghysels ierr = VecRestoreArray(x,&xptr);CHKERRQ(ierr); 2337d6ea485SPieter Ghysels ierr = VecRestoreArrayRead(b_mpi,&bptr);CHKERRQ(ierr); 2347d6ea485SPieter Ghysels PetscFunctionReturn(0); 2357d6ea485SPieter Ghysels } 2367d6ea485SPieter Ghysels 237*291d0ed5SPieter Ghysels #undef __FUNCT__ 238*291d0ed5SPieter Ghysels #define __FUNCT__ "MatMatSolve_STRUMPACK" 239ad0c5e61SPieter Ghysels static PetscErrorCode MatMatSolve_STRUMPACK(Mat A,Mat B_mpi,Mat X) 2407d6ea485SPieter Ghysels { 2417d6ea485SPieter Ghysels PetscErrorCode ierr; 2427d6ea485SPieter Ghysels PetscBool flg; 2437d6ea485SPieter Ghysels 2447d6ea485SPieter Ghysels PetscFunctionBegin; 2457d6ea485SPieter Ghysels ierr = PetscObjectTypeCompareAny((PetscObject)B_mpi,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 2467d6ea485SPieter Ghysels if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 2477d6ea485SPieter Ghysels ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 2487d6ea485SPieter Ghysels if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); 2497d6ea485SPieter Ghysels SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_STRUMPACK() is not implemented yet"); 2507d6ea485SPieter Ghysels PetscFunctionReturn(0); 2517d6ea485SPieter Ghysels } 2527d6ea485SPieter Ghysels 253*291d0ed5SPieter Ghysels #undef __FUNCT__ 254*291d0ed5SPieter Ghysels #define __FUNCT__ "MatFactorInfo_STRUMPACK" 255ad0c5e61SPieter Ghysels static PetscErrorCode MatFactorInfo_STRUMPACK(Mat A,PetscViewer viewer) 256ad0c5e61SPieter Ghysels { 257ad0c5e61SPieter Ghysels PetscErrorCode ierr; 258ad0c5e61SPieter Ghysels 259ad0c5e61SPieter Ghysels PetscFunctionBegin; 260ad0c5e61SPieter Ghysels /* check if matrix is strumpack type */ 261ad0c5e61SPieter Ghysels if (A->ops->solve != MatSolve_STRUMPACK) PetscFunctionReturn(0); 262ad0c5e61SPieter Ghysels ierr = PetscViewerASCIIPrintf(viewer,"STRUMPACK sparse solver!\n");CHKERRQ(ierr); 263ad0c5e61SPieter Ghysels PetscFunctionReturn(0); 264ad0c5e61SPieter Ghysels } 265ad0c5e61SPieter Ghysels 266*291d0ed5SPieter Ghysels #undef __FUNCT__ 267*291d0ed5SPieter Ghysels #define __FUNCT__ "MatView_STRUMPACK" 268ad0c5e61SPieter Ghysels static PetscErrorCode MatView_STRUMPACK(Mat A,PetscViewer viewer) 269ad0c5e61SPieter Ghysels { 270ad0c5e61SPieter Ghysels PetscErrorCode ierr; 271ad0c5e61SPieter Ghysels PetscBool iascii; 272ad0c5e61SPieter Ghysels PetscViewerFormat format; 273ad0c5e61SPieter Ghysels 274ad0c5e61SPieter Ghysels PetscFunctionBegin; 275ad0c5e61SPieter Ghysels ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 276ad0c5e61SPieter Ghysels if (iascii) { 277ad0c5e61SPieter Ghysels ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 278ad0c5e61SPieter Ghysels if (format == PETSC_VIEWER_ASCII_INFO) { 279ad0c5e61SPieter Ghysels ierr = MatFactorInfo_STRUMPACK(A,viewer);CHKERRQ(ierr); 280ad0c5e61SPieter Ghysels } 281ad0c5e61SPieter Ghysels } 282ad0c5e61SPieter Ghysels PetscFunctionReturn(0); 283ad0c5e61SPieter Ghysels } 2847d6ea485SPieter Ghysels 285*291d0ed5SPieter Ghysels #undef __FUNCT__ 286*291d0ed5SPieter Ghysels #define __FUNCT__ "MatLUFactorNumeric_STRUMPACK" 287ad0c5e61SPieter Ghysels static PetscErrorCode MatLUFactorNumeric_STRUMPACK(Mat F,Mat A,const MatFactorInfo *info) 2887d6ea485SPieter Ghysels { 2890d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 2907d6ea485SPieter Ghysels STRUMPACK_RETURN_CODE sp_err; 2917d6ea485SPieter Ghysels Mat_SeqAIJ *A_d,*A_o; 2927d6ea485SPieter Ghysels Mat_MPIAIJ *mat; 2937d6ea485SPieter Ghysels PetscErrorCode ierr; 2940d08a34dSPieter Ghysels PetscInt M=A->rmap->N,m=A->rmap->n; 2957d6ea485SPieter Ghysels PetscBool flg; 2967d6ea485SPieter Ghysels 2977d6ea485SPieter Ghysels PetscFunctionBegin; 2987d6ea485SPieter Ghysels ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); 2997d6ea485SPieter Ghysels if (flg) { /* A is MATMPIAIJ */ 3007d6ea485SPieter Ghysels mat = (Mat_MPIAIJ*)A->data; 3017d6ea485SPieter Ghysels A_d = (Mat_SeqAIJ*)(mat->A)->data; 3027d6ea485SPieter Ghysels A_o = (Mat_SeqAIJ*)(mat->B)->data; 3030d08a34dSPieter 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)); 3047d6ea485SPieter Ghysels } else { /* A is MATSEQAIJ */ 3057d6ea485SPieter Ghysels A_d = (Mat_SeqAIJ*)A->data; 3060d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_set_csr_matrix",STRUMPACK_set_csr_matrix(*S,&M,A_d->i,A_d->j,A_d->a,0)); 3077d6ea485SPieter Ghysels } 3087d6ea485SPieter Ghysels 3097d6ea485SPieter Ghysels /* Reorder and Factor the matrix. */ 3107d6ea485SPieter Ghysels /* TODO figure out how to avoid reorder if the matrix values changed, but the pattern remains the same. */ 3110d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_reorder",sp_err = STRUMPACK_reorder(*S)); 3120d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_factor",sp_err = STRUMPACK_factor(*S)); 3130d08a34dSPieter Ghysels switch (sp_err) { 3140d08a34dSPieter Ghysels case STRUMPACK_SUCCESS: break; 3150d08a34dSPieter Ghysels case STRUMPACK_MATRIX_NOT_SET: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); break; } 3160d08a34dSPieter Ghysels case STRUMPACK_REORDERING_ERROR: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); break; } 3170d08a34dSPieter Ghysels default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: factorization failed"); 3187d6ea485SPieter Ghysels } 3197d6ea485SPieter Ghysels PetscFunctionReturn(0); 3207d6ea485SPieter Ghysels } 3217d6ea485SPieter Ghysels 322*291d0ed5SPieter Ghysels #undef __FUNCT__ 323*291d0ed5SPieter Ghysels #define __FUNCT__ "MatLUFactorSymbolic_STRUMPACK" 324ad0c5e61SPieter Ghysels static PetscErrorCode MatLUFactorSymbolic_STRUMPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 3257d6ea485SPieter Ghysels { 3267d6ea485SPieter Ghysels PetscFunctionBegin; 3277d6ea485SPieter Ghysels F->ops->lufactornumeric = MatLUFactorNumeric_STRUMPACK; 3287d6ea485SPieter Ghysels F->ops->solve = MatSolve_STRUMPACK; 3297d6ea485SPieter Ghysels F->ops->matsolve = MatMatSolve_STRUMPACK; 3307d6ea485SPieter Ghysels PetscFunctionReturn(0); 3317d6ea485SPieter Ghysels } 3327d6ea485SPieter Ghysels 333*291d0ed5SPieter Ghysels #undef __FUNCT__ 334*291d0ed5SPieter Ghysels #define __FUNCT__ "MatFactorGetSolverPackage_aij_strumpack" 335ad0c5e61SPieter Ghysels static PetscErrorCode MatFactorGetSolverPackage_aij_strumpack(Mat A,const MatSolverPackage *type) 3367d6ea485SPieter Ghysels { 3377d6ea485SPieter Ghysels PetscFunctionBegin; 3387d6ea485SPieter Ghysels *type = MATSOLVERSTRUMPACK; 3397d6ea485SPieter Ghysels PetscFunctionReturn(0); 3407d6ea485SPieter Ghysels } 3417d6ea485SPieter Ghysels 342575704cbSPieter Ghysels /*MC 343575704cbSPieter Ghysels MATSOLVERSSTRUMPACK = "strumpack" - A solver package providing a direct sparse solver (PCLU) 344575704cbSPieter Ghysels and a preconditioner (PCILU) using low-rank compression via the external package STRUMPACK. 345575704cbSPieter Ghysels 346575704cbSPieter Ghysels Consult the STRUMPACK-sparse manual for more info. 347575704cbSPieter Ghysels 348575704cbSPieter Ghysels Use 349575704cbSPieter Ghysels ./configure --download-strumpack 350575704cbSPieter Ghysels to have PETSc installed with STRUMPACK 351575704cbSPieter Ghysels 352575704cbSPieter Ghysels Use 353575704cbSPieter Ghysels -pc_type lu -pc_factor_mat_solver_package strumpack 354575704cbSPieter Ghysels to use this as an exact (direct) solver, use 355575704cbSPieter Ghysels -pc_type ilu -pc_factor_mat_solver_package strumpack 356575704cbSPieter Ghysels to enable low-rank compression (i.e, use as a preconditioner). 357575704cbSPieter Ghysels 358575704cbSPieter Ghysels Works with AIJ matrices 359575704cbSPieter Ghysels 360575704cbSPieter Ghysels Options Database Keys: 361575704cbSPieter Ghysels + -mat_strumpack_rctol <1e-4> - Relative compression tolerance (None) 362*291d0ed5SPieter Ghysels . -mat_strumpack_hss_min_front_size - Minimum size of dense block for HSS compression (None) 363*291d0ed5SPieter Ghysels . -mat_strumpack_hss_min_sep_size - Minimum size of separator for HSS compression (None) 364575704cbSPieter Ghysels . -mat_strumpack_colperm <TRUE> - Permute matrix to make diagonal nonzeros (None) 365575704cbSPieter Ghysels 366575704cbSPieter Ghysels Level: beginner 367575704cbSPieter Ghysels 368575704cbSPieter Ghysels .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, PCFactorSetMatSolverPackage(), MatSolverPackage 369575704cbSPieter Ghysels M*/ 370*291d0ed5SPieter Ghysels #undef __FUNCT__ 371*291d0ed5SPieter Ghysels #define __FUNCT__ "MatGetFactor_aij_strumpack" 372ad0c5e61SPieter Ghysels static PetscErrorCode MatGetFactor_aij_strumpack(Mat A,MatFactorType ftype,Mat *F) 3737d6ea485SPieter Ghysels { 3747d6ea485SPieter Ghysels Mat B; 3750d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S; 3767d6ea485SPieter Ghysels PetscErrorCode ierr; 3777d6ea485SPieter Ghysels PetscInt M=A->rmap->N,N=A->cmap->N; 3787d6ea485SPieter Ghysels STRUMPACK_INTERFACE iface; 379575704cbSPieter Ghysels PetscBool verb,flg,set; 380575704cbSPieter Ghysels PetscReal rctol; 381575704cbSPieter Ghysels PetscInt hssminsize; 3827d6ea485SPieter Ghysels int argc; 38355c022e5SPieter Ghysels char **args,*copts,*pname; 38455c022e5SPieter Ghysels size_t len; 3857d6ea485SPieter Ghysels const STRUMPACK_PRECISION table[2][2][2] = {{{STRUMPACK_FLOATCOMPLEX_64,STRUMPACK_DOUBLECOMPLEX_64}, 3867d6ea485SPieter Ghysels {STRUMPACK_FLOAT_64,STRUMPACK_DOUBLE_64}}, 3877d6ea485SPieter Ghysels {{STRUMPACK_FLOATCOMPLEX,STRUMPACK_DOUBLECOMPLEX}, 3887d6ea485SPieter Ghysels {STRUMPACK_FLOAT,STRUMPACK_DOUBLE}}}; 3897d6ea485SPieter Ghysels const STRUMPACK_PRECISION prec = table[(sizeof(PetscInt)==8)?0:1][(PETSC_SCALAR==PETSC_COMPLEX)?0:1][(PETSC_REAL==PETSC_FLOAT)?0:1]; 3907d6ea485SPieter Ghysels 3917d6ea485SPieter Ghysels PetscFunctionBegin; 3927d6ea485SPieter Ghysels /* Create the factorization matrix */ 3937d6ea485SPieter Ghysels ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 3947d6ea485SPieter Ghysels ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,M,N);CHKERRQ(ierr); 3957d6ea485SPieter Ghysels ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 3967d6ea485SPieter Ghysels ierr = MatSeqAIJSetPreallocation(B,0,NULL); 3977d6ea485SPieter Ghysels ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr); 398575704cbSPieter Ghysels if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) { 3997d6ea485SPieter Ghysels B->ops->lufactorsymbolic = MatLUFactorSymbolic_STRUMPACK; 400575704cbSPieter Ghysels B->ops->ilufactorsymbolic = MatLUFactorSymbolic_STRUMPACK; 401575704cbSPieter Ghysels } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 4027d6ea485SPieter Ghysels B->ops->view = MatView_STRUMPACK; 4037d6ea485SPieter Ghysels B->ops->destroy = MatDestroy_STRUMPACK; 4047d6ea485SPieter Ghysels B->ops->getdiagonal = MatGetDiagonal_STRUMPACK; 4057d6ea485SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_aij_strumpack);CHKERRQ(ierr); 406575704cbSPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSRelCompTol_C",MatSTRUMPACKSetHSSRelCompTol_STRUMPACK);CHKERRQ(ierr); 407*291d0ed5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSMinFrontSize_C",MatSTRUMPACKSetHSSMinFrontSize_STRUMPACK);CHKERRQ(ierr); 408*291d0ed5SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSMinSepSize_C",MatSTRUMPACKSetHSSMinSepSize_STRUMPACK);CHKERRQ(ierr); 409575704cbSPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetColPerm_C",MatSTRUMPACKSetColPerm_STRUMPACK);CHKERRQ(ierr); 410575704cbSPieter Ghysels B->factortype = ftype; 4110d08a34dSPieter Ghysels ierr = PetscNewLog(B,&S);CHKERRQ(ierr); 4120d08a34dSPieter Ghysels B->spptr = S; 4130d08a34dSPieter Ghysels 4140d08a34dSPieter Ghysels ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg); 4150d08a34dSPieter Ghysels iface = flg ? STRUMPACK_MT : STRUMPACK_MPI_DIST; 4167d6ea485SPieter Ghysels 4177d6ea485SPieter Ghysels ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"STRUMPACK Options","Mat");CHKERRQ(ierr); 4187d6ea485SPieter Ghysels 419575704cbSPieter Ghysels verb = PetscLogPrintInfo ? PETSC_TRUE : PETSC_FALSE; 4207d6ea485SPieter Ghysels ierr = PetscOptionsBool("-mat_strumpack_verbose","Print STRUMPACK information","None",verb,&verb,NULL);CHKERRQ(ierr); 4217d6ea485SPieter Ghysels 42255c022e5SPieter Ghysels ierr = PetscOptionsGetAll(NULL,&copts);CHKERRQ(ierr); 42355c022e5SPieter Ghysels ierr = PetscStrlen(copts,&len);CHKERRQ(ierr); 42455c022e5SPieter Ghysels len += PETSC_MAX_PATH_LEN+1; 42555c022e5SPieter Ghysels ierr = PetscMalloc1(len,&pname);CHKERRQ(ierr); 42655c022e5SPieter Ghysels /* first string is assumed to be the program name, so add program name to options string */ 42755c022e5SPieter Ghysels ierr = PetscGetProgramName(pname,len);CHKERRQ(ierr); 42855c022e5SPieter Ghysels ierr = PetscStrcat(pname," ");CHKERRQ(ierr); 42955c022e5SPieter Ghysels ierr = PetscStrcat(pname,copts);CHKERRQ(ierr); 430ee5a0f3aSPieter Ghysels ierr = PetscStrToArray(pname,' ',&argc,&args);CHKERRQ(ierr); 43155c022e5SPieter Ghysels ierr = PetscFree(copts);CHKERRQ(ierr); 43255c022e5SPieter Ghysels ierr = PetscFree(pname);CHKERRQ(ierr); 43355c022e5SPieter Ghysels 4340d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_init",STRUMPACK_init(S,PetscObjectComm((PetscObject)A),prec,iface,argc,args,verb)); 4357d6ea485SPieter Ghysels 436*291d0ed5SPieter Ghysels PetscStackCall("STRUMPACK_hss_rel_tol",rctol = (PetscReal)STRUMPACK_hss_rel_tol(*S)); 437*291d0ed5SPieter Ghysels ierr = PetscOptionsReal("-mat_strumpack_hss_rel_tol","Relative compression tolerance","None",rctol,&rctol,&set);CHKERRQ(ierr); 438*291d0ed5SPieter Ghysels if (set) PetscStackCall("STRUMPACK_set_hss_rel_tol",STRUMPACK_set_hss_rel_tol(*S,(double)rctol)); 439575704cbSPieter Ghysels 440*291d0ed5SPieter Ghysels PetscStackCall("STRUMPACK_mc64job",flg = (STRUMPACK_mc64job(*S) == 0) ? PETSC_FALSE : PETSC_TRUE); 441575704cbSPieter Ghysels ierr = PetscOptionsBool("-mat_strumpack_colperm","Find a col perm to get nonzero diagonal","None",flg,&flg,&set);CHKERRQ(ierr); 4420d08a34dSPieter Ghysels if (set) PetscStackCall("STRUMPACK_set_mc64job",STRUMPACK_set_mc64job(*S,flg ? 5 : 0)); 443575704cbSPieter Ghysels 444*291d0ed5SPieter Ghysels PetscStackCall("STRUMPACK_HSS_min_front_size",hssminsize = (PetscInt)STRUMPACK_HSS_min_front_size(*S)); 445*291d0ed5SPieter Ghysels ierr = PetscOptionsInt("-mat_strumpack_hss_min_front_size","Minimum size of dense block for HSS compression","None",hssminsize,&hssminsize,&set);CHKERRQ(ierr); 446*291d0ed5SPieter Ghysels if (set) PetscStackCall("STRUMPACK_set_HSS_min_front_size",STRUMPACK_set_HSS_min_front_size(*S,(int)hssminsize)); 447*291d0ed5SPieter Ghysels 448*291d0ed5SPieter Ghysels PetscStackCall("STRUMPACK_HSS_min_sep_size",hssminsize = (PetscInt)STRUMPACK_HSS_min_sep_size(*S)); 449*291d0ed5SPieter Ghysels ierr = PetscOptionsInt("-mat_strumpack_hss_min_sep_size","Minimum size of separator for HSS compression","None",hssminsize,&hssminsize,&set);CHKERRQ(ierr); 450*291d0ed5SPieter Ghysels if (set) PetscStackCall("STRUMPACK_set_HSS_min_sep_size",STRUMPACK_set_HSS_min_sep_size(*S,(int)hssminsize)); 451575704cbSPieter Ghysels 452575704cbSPieter Ghysels PetscOptionsEnd(); 453575704cbSPieter Ghysels 45488382b05SPieter Ghysels if (ftype == MAT_FACTOR_ILU) { 45588382b05SPieter Ghysels /* When enabling HSS compression, the STRUMPACK solver becomes an incomplete */ 45688382b05SPieter Ghysels /* (or approximate) LU factorization. */ 457*291d0ed5SPieter Ghysels PetscStackCall("STRUMPACK_enable_HSS",STRUMPACK_enable_HSS(*S)); 458575704cbSPieter Ghysels /* Disable the outer iterative solver from STRUMPACK. */ 459575704cbSPieter Ghysels /* When STRUMPACK is used as a direct solver, it will by default do iterative refinement. */ 460575704cbSPieter Ghysels /* When STRUMPACK is used with as an approximate factorization preconditioner (by enabling */ 46188382b05SPieter Ghysels /* low-rank compression), it will use it's own GMRES. Here we can disable the */ 462575704cbSPieter Ghysels /* outer iterative solver, as PETSc uses STRUMPACK from within a KSP. */ 4630d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_set_Krylov_solver", STRUMPACK_set_Krylov_solver(*S, STRUMPACK_DIRECT)); 46488382b05SPieter Ghysels } 465575704cbSPieter Ghysels 466575704cbSPieter Ghysels /* Allow the user to set or overwrite the above options from the command line */ 4670d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_set_from_options",STRUMPACK_set_from_options(*S)); 46855c022e5SPieter Ghysels ierr = PetscStrToArrayDestroy(argc,args);CHKERRQ(ierr); 46955c022e5SPieter Ghysels 4707d6ea485SPieter Ghysels *F = B; 4717d6ea485SPieter Ghysels PetscFunctionReturn(0); 4727d6ea485SPieter Ghysels } 4737d6ea485SPieter Ghysels 474*291d0ed5SPieter Ghysels #undef __FUNCT__ 475*291d0ed5SPieter Ghysels #define __FUNCT__ "MatSolverPackageRegister_STRUMPACK" 4767d6ea485SPieter Ghysels PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_STRUMPACK(void) 4777d6ea485SPieter Ghysels { 4787d6ea485SPieter Ghysels PetscErrorCode ierr; 4797d6ea485SPieter Ghysels 4807d6ea485SPieter Ghysels PetscFunctionBegin; 4817d6ea485SPieter Ghysels ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr); 4827d6ea485SPieter Ghysels ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr); 483575704cbSPieter Ghysels ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_ILU,MatGetFactor_aij_strumpack);CHKERRQ(ierr); 484575704cbSPieter Ghysels ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_ILU,MatGetFactor_aij_strumpack);CHKERRQ(ierr); 4857d6ea485SPieter Ghysels PetscFunctionReturn(0); 4867d6ea485SPieter Ghysels } 487