17d6ea485SPieter Ghysels #include <../src/mat/impls/aij/seq/aij.h> 27d6ea485SPieter Ghysels #include <../src/mat/impls/aij/mpi/mpiaij.h> 37d6ea485SPieter Ghysels #include <StrumpackSparseSolver.h> 47d6ea485SPieter Ghysels 57d6ea485SPieter Ghysels #undef __FUNCT__ 67d6ea485SPieter 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 147d6ea485SPieter Ghysels #undef __FUNCT__ 157d6ea485SPieter Ghysels #define __FUNCT__ "MatDestroy_STRUMPACK" 16ad0c5e61SPieter Ghysels static PetscErrorCode MatDestroy_STRUMPACK(Mat A) 177d6ea485SPieter Ghysels { 18*0d08a34dSPieter 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 */ 24*0d08a34dSPieter 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); 36575704cbSPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSMinSize_C",NULL);CHKERRQ(ierr); 37575704cbSPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetColPerm_C",NULL);CHKERRQ(ierr); 38575704cbSPieter Ghysels 39575704cbSPieter Ghysels PetscFunctionReturn(0); 40575704cbSPieter Ghysels } 41575704cbSPieter Ghysels 42575704cbSPieter Ghysels #undef __FUNCT__ 43575704cbSPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSRelCompTol_STRUMPACK" 44575704cbSPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSRelCompTol_STRUMPACK(Mat F,PetscReal rctol) 45575704cbSPieter Ghysels { 46*0d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 47575704cbSPieter Ghysels 48575704cbSPieter Ghysels PetscFunctionBegin; 49*0d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_set_rctol", STRUMPACK_set_rctol(*S,rctol)); 50575704cbSPieter Ghysels PetscFunctionReturn(0); 51575704cbSPieter Ghysels } 52575704cbSPieter Ghysels 53575704cbSPieter Ghysels #undef __FUNCT__ 54575704cbSPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSRelCompTol" 55575704cbSPieter Ghysels /*@ 56575704cbSPieter Ghysels MatSTRUMPACKSetHSSRelCompTol - Set STRUMPACK relative tolerance for HSS compression 57575704cbSPieter Ghysels Logically Collective on Mat 58575704cbSPieter Ghysels 59575704cbSPieter Ghysels Input Parameters: 60575704cbSPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 61575704cbSPieter Ghysels - rctol - relative compression tolerance 62575704cbSPieter Ghysels 63575704cbSPieter Ghysels Options Database: 64575704cbSPieter Ghysels . -mat_strumpack_rctol <rctol> 65575704cbSPieter Ghysels 66575704cbSPieter Ghysels Level: beginner 67575704cbSPieter Ghysels 68575704cbSPieter Ghysels References: 69575704cbSPieter Ghysels . STRUMPACK manual 70575704cbSPieter Ghysels 71575704cbSPieter Ghysels .seealso: MatGetFactor() 72575704cbSPieter Ghysels @*/ 73575704cbSPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat F,PetscReal rctol) 74575704cbSPieter Ghysels { 75575704cbSPieter Ghysels PetscErrorCode ierr; 76575704cbSPieter Ghysels 77575704cbSPieter Ghysels PetscFunctionBegin; 78575704cbSPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 79575704cbSPieter Ghysels PetscValidLogicalCollectiveReal(F,rctol,2); 80575704cbSPieter Ghysels ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSRelCompTol_C",(Mat,PetscReal),(F,rctol));CHKERRQ(ierr); 81575704cbSPieter Ghysels PetscFunctionReturn(0); 82575704cbSPieter Ghysels } 83575704cbSPieter Ghysels 84575704cbSPieter Ghysels #undef __FUNCT__ 85575704cbSPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSMinSize_STRUMPACK" 86575704cbSPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSMinSize_STRUMPACK(Mat F,PetscInt hssminsize) 87575704cbSPieter Ghysels { 88*0d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 89575704cbSPieter Ghysels 90575704cbSPieter Ghysels PetscFunctionBegin; 91*0d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_set_minimum_HSS_size", STRUMPACK_set_minimum_HSS_size(*S,hssminsize)); 92575704cbSPieter Ghysels PetscFunctionReturn(0); 93575704cbSPieter Ghysels } 94575704cbSPieter Ghysels 95575704cbSPieter Ghysels #undef __FUNCT__ 96575704cbSPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSMinSize" 97575704cbSPieter Ghysels /*@ 98575704cbSPieter Ghysels MatSTRUMPACKSetHSSMinSize - Set STRUMPACK minimum dense matrix size for low-rank approximation 99575704cbSPieter Ghysels Logically Collective on Mat 100575704cbSPieter Ghysels 101575704cbSPieter Ghysels Input Parameters: 102575704cbSPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 103575704cbSPieter Ghysels - hssminsize - minimum dense matrix size for low-rank approximation 104575704cbSPieter Ghysels 105575704cbSPieter Ghysels Options Database: 106575704cbSPieter Ghysels . -mat_strumpack_hssminsize <hssminsize> 107575704cbSPieter Ghysels 108575704cbSPieter Ghysels Level: beginner 109575704cbSPieter Ghysels 110575704cbSPieter Ghysels References: 111575704cbSPieter Ghysels . STRUMPACK manual 112575704cbSPieter Ghysels 113575704cbSPieter Ghysels .seealso: MatGetFactor() 114575704cbSPieter Ghysels @*/ 115575704cbSPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSMinSize(Mat F,PetscInt hssminsize) 116575704cbSPieter Ghysels { 117575704cbSPieter Ghysels PetscErrorCode ierr; 118575704cbSPieter Ghysels 119575704cbSPieter Ghysels PetscFunctionBegin; 120575704cbSPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 121575704cbSPieter Ghysels PetscValidLogicalCollectiveInt(F,hssminsize,2); 122575704cbSPieter Ghysels ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSMinSize_C",(Mat,PetscInt),(F,hssminsize));CHKERRQ(ierr); 123575704cbSPieter Ghysels PetscFunctionReturn(0); 124575704cbSPieter Ghysels } 125575704cbSPieter Ghysels 126575704cbSPieter Ghysels #undef __FUNCT__ 127575704cbSPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetColPerm_STRUMPACK" 128575704cbSPieter Ghysels static PetscErrorCode MatSTRUMPACKSetColPerm_STRUMPACK(Mat F,PetscBool cperm) 129575704cbSPieter Ghysels { 130*0d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 131575704cbSPieter Ghysels 132575704cbSPieter Ghysels PetscFunctionBegin; 133*0d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_set_mc64job",STRUMPACK_set_mc64job(*S,cperm ? 5 : 0)); 134575704cbSPieter Ghysels PetscFunctionReturn(0); 135575704cbSPieter Ghysels } 136575704cbSPieter Ghysels 137575704cbSPieter Ghysels #undef __FUNCT__ 138575704cbSPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetColPerm" 139575704cbSPieter Ghysels /*@ 140575704cbSPieter Ghysels MatSTRUMPACKSetColPerm - Set whether STRUMPACK should try to permute the columns of the matrix in order to get a nonzero diagonal 141575704cbSPieter Ghysels Logically Collective on Mat 142575704cbSPieter Ghysels 143575704cbSPieter Ghysels Input Parameters: 144575704cbSPieter Ghysels + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 145575704cbSPieter Ghysels - cperm - whether or not to permute (internally) the columns of the matrix 146575704cbSPieter Ghysels 147575704cbSPieter Ghysels Options Database: 148575704cbSPieter Ghysels . -mat_strumpack_colperm <cperm> 149575704cbSPieter Ghysels 150575704cbSPieter Ghysels Level: beginner 151575704cbSPieter Ghysels 152575704cbSPieter Ghysels References: 153575704cbSPieter Ghysels . STRUMPACK manual 154575704cbSPieter Ghysels 155575704cbSPieter Ghysels .seealso: MatGetFactor() 156575704cbSPieter Ghysels @*/ 157575704cbSPieter Ghysels PetscErrorCode MatSTRUMPACKSetColPerm(Mat F,PetscBool cperm) 158575704cbSPieter Ghysels { 159575704cbSPieter Ghysels PetscErrorCode ierr; 160575704cbSPieter Ghysels 161575704cbSPieter Ghysels PetscFunctionBegin; 162575704cbSPieter Ghysels PetscValidHeaderSpecific(F,MAT_CLASSID,1); 163575704cbSPieter Ghysels PetscValidLogicalCollectiveBool(F,cperm,2); 164575704cbSPieter Ghysels ierr = PetscTryMethod(F,"MatSTRUMPACKSetColPerm_C",(Mat,PetscBool),(F,cperm));CHKERRQ(ierr); 1657d6ea485SPieter Ghysels PetscFunctionReturn(0); 1667d6ea485SPieter Ghysels } 1677d6ea485SPieter Ghysels 1687d6ea485SPieter Ghysels #undef __FUNCT__ 1697d6ea485SPieter Ghysels #define __FUNCT__ "MatSolve_STRUMPACK" 170ad0c5e61SPieter Ghysels static PetscErrorCode MatSolve_STRUMPACK(Mat A,Vec b_mpi,Vec x) 1717d6ea485SPieter Ghysels { 172*0d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)A->spptr; 1737d6ea485SPieter Ghysels STRUMPACK_RETURN_CODE sp_err; 1747d6ea485SPieter Ghysels PetscErrorCode ierr; 1757d6ea485SPieter Ghysels const PetscScalar *bptr; 1767d6ea485SPieter Ghysels PetscScalar *xptr; 1777d6ea485SPieter Ghysels 1787d6ea485SPieter Ghysels PetscFunctionBegin; 1797d6ea485SPieter Ghysels ierr = VecGetArray(x,&xptr);CHKERRQ(ierr); 1807d6ea485SPieter Ghysels ierr = VecGetArrayRead(b_mpi,&bptr);CHKERRQ(ierr); 181*0d08a34dSPieter Ghysels 182*0d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_solve",sp_err = STRUMPACK_solve(*S,(PetscScalar*)bptr,xptr,0)); 183*0d08a34dSPieter Ghysels switch (sp_err) { 184*0d08a34dSPieter Ghysels case STRUMPACK_SUCCESS: break; 185*0d08a34dSPieter Ghysels case STRUMPACK_MATRIX_NOT_SET: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); break; } 186*0d08a34dSPieter Ghysels case STRUMPACK_REORDERING_ERROR: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); break; } 187*0d08a34dSPieter Ghysels default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: solve failed"); 1887d6ea485SPieter Ghysels } 1897d6ea485SPieter Ghysels ierr = VecRestoreArray(x,&xptr);CHKERRQ(ierr); 1907d6ea485SPieter Ghysels ierr = VecRestoreArrayRead(b_mpi,&bptr);CHKERRQ(ierr); 1917d6ea485SPieter Ghysels PetscFunctionReturn(0); 1927d6ea485SPieter Ghysels } 1937d6ea485SPieter Ghysels 1947d6ea485SPieter Ghysels #undef __FUNCT__ 1957d6ea485SPieter Ghysels #define __FUNCT__ "MatMatSolve_STRUMPACK" 196ad0c5e61SPieter Ghysels static PetscErrorCode MatMatSolve_STRUMPACK(Mat A,Mat B_mpi,Mat X) 1977d6ea485SPieter Ghysels { 1987d6ea485SPieter Ghysels PetscErrorCode ierr; 1997d6ea485SPieter Ghysels PetscBool flg; 2007d6ea485SPieter Ghysels 2017d6ea485SPieter Ghysels PetscFunctionBegin; 2027d6ea485SPieter Ghysels ierr = PetscObjectTypeCompareAny((PetscObject)B_mpi,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 2037d6ea485SPieter Ghysels if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 2047d6ea485SPieter Ghysels ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 2057d6ea485SPieter Ghysels if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); 2067d6ea485SPieter Ghysels SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_STRUMPACK() is not implemented yet"); 2077d6ea485SPieter Ghysels PetscFunctionReturn(0); 2087d6ea485SPieter Ghysels } 2097d6ea485SPieter Ghysels 210ad0c5e61SPieter Ghysels #undef __FUNCT__ 211ad0c5e61SPieter Ghysels #define __FUNCT__ "MatFactorInfo_STRUMPACK" 212ad0c5e61SPieter Ghysels static PetscErrorCode MatFactorInfo_STRUMPACK(Mat A,PetscViewer viewer) 213ad0c5e61SPieter Ghysels { 214ad0c5e61SPieter Ghysels PetscErrorCode ierr; 215ad0c5e61SPieter Ghysels 216ad0c5e61SPieter Ghysels PetscFunctionBegin; 217ad0c5e61SPieter Ghysels /* check if matrix is strumpack type */ 218ad0c5e61SPieter Ghysels if (A->ops->solve != MatSolve_STRUMPACK) PetscFunctionReturn(0); 219ad0c5e61SPieter Ghysels ierr = PetscViewerASCIIPrintf(viewer,"STRUMPACK sparse solver!\n");CHKERRQ(ierr); 220ad0c5e61SPieter Ghysels PetscFunctionReturn(0); 221ad0c5e61SPieter Ghysels } 222ad0c5e61SPieter Ghysels 223ad0c5e61SPieter Ghysels #undef __FUNCT__ 224ad0c5e61SPieter Ghysels #define __FUNCT__ "MatView_STRUMPACK" 225ad0c5e61SPieter Ghysels static PetscErrorCode MatView_STRUMPACK(Mat A,PetscViewer viewer) 226ad0c5e61SPieter Ghysels { 227ad0c5e61SPieter Ghysels PetscErrorCode ierr; 228ad0c5e61SPieter Ghysels PetscBool iascii; 229ad0c5e61SPieter Ghysels PetscViewerFormat format; 230ad0c5e61SPieter Ghysels 231ad0c5e61SPieter Ghysels PetscFunctionBegin; 232ad0c5e61SPieter Ghysels ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 233ad0c5e61SPieter Ghysels if (iascii) { 234ad0c5e61SPieter Ghysels ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 235ad0c5e61SPieter Ghysels if (format == PETSC_VIEWER_ASCII_INFO) { 236ad0c5e61SPieter Ghysels ierr = MatFactorInfo_STRUMPACK(A,viewer);CHKERRQ(ierr); 237ad0c5e61SPieter Ghysels } 238ad0c5e61SPieter Ghysels } 239ad0c5e61SPieter Ghysels PetscFunctionReturn(0); 240ad0c5e61SPieter Ghysels } 2417d6ea485SPieter Ghysels 2427d6ea485SPieter Ghysels #undef __FUNCT__ 2437d6ea485SPieter Ghysels #define __FUNCT__ "MatLUFactorNumeric_STRUMPACK" 244ad0c5e61SPieter Ghysels static PetscErrorCode MatLUFactorNumeric_STRUMPACK(Mat F,Mat A,const MatFactorInfo *info) 2457d6ea485SPieter Ghysels { 246*0d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 2477d6ea485SPieter Ghysels STRUMPACK_RETURN_CODE sp_err; 2487d6ea485SPieter Ghysels Mat_SeqAIJ *A_d,*A_o; 2497d6ea485SPieter Ghysels Mat_MPIAIJ *mat; 2507d6ea485SPieter Ghysels PetscErrorCode ierr; 251*0d08a34dSPieter Ghysels PetscInt M=A->rmap->N,m=A->rmap->n; 2527d6ea485SPieter Ghysels PetscBool flg; 2537d6ea485SPieter Ghysels 2547d6ea485SPieter Ghysels PetscFunctionBegin; 2557d6ea485SPieter Ghysels ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); 2567d6ea485SPieter Ghysels if (flg) { /* A is MATMPIAIJ */ 2577d6ea485SPieter Ghysels mat = (Mat_MPIAIJ*)A->data; 2587d6ea485SPieter Ghysels A_d = (Mat_SeqAIJ*)(mat->A)->data; 2597d6ea485SPieter Ghysels A_o = (Mat_SeqAIJ*)(mat->B)->data; 260*0d08a34dSPieter 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)); 2617d6ea485SPieter Ghysels } else { /* A is MATSEQAIJ */ 2627d6ea485SPieter Ghysels A_d = (Mat_SeqAIJ*)A->data; 263*0d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_set_csr_matrix",STRUMPACK_set_csr_matrix(*S,&M,A_d->i,A_d->j,A_d->a,0)); 2647d6ea485SPieter Ghysels } 2657d6ea485SPieter Ghysels 2667d6ea485SPieter Ghysels /* Reorder and Factor the matrix. */ 2677d6ea485SPieter Ghysels /* TODO figure out how to avoid reorder if the matrix values changed, but the pattern remains the same. */ 268*0d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_reorder",sp_err = STRUMPACK_reorder(*S)); 269*0d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_factor",sp_err = STRUMPACK_factor(*S)); 270*0d08a34dSPieter Ghysels switch (sp_err) { 271*0d08a34dSPieter Ghysels case STRUMPACK_SUCCESS: break; 272*0d08a34dSPieter Ghysels case STRUMPACK_MATRIX_NOT_SET: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); break; } 273*0d08a34dSPieter Ghysels case STRUMPACK_REORDERING_ERROR: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); break; } 274*0d08a34dSPieter Ghysels default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: factorization failed"); 2757d6ea485SPieter Ghysels } 2767d6ea485SPieter Ghysels PetscFunctionReturn(0); 2777d6ea485SPieter Ghysels } 2787d6ea485SPieter Ghysels 2797d6ea485SPieter Ghysels #undef __FUNCT__ 2807d6ea485SPieter Ghysels #define __FUNCT__ "MatLUFactorSymbolic_STRUMPACK" 281ad0c5e61SPieter Ghysels static PetscErrorCode MatLUFactorSymbolic_STRUMPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 2827d6ea485SPieter Ghysels { 2837d6ea485SPieter Ghysels PetscFunctionBegin; 2847d6ea485SPieter Ghysels F->ops->lufactornumeric = MatLUFactorNumeric_STRUMPACK; 2857d6ea485SPieter Ghysels F->ops->solve = MatSolve_STRUMPACK; 2867d6ea485SPieter Ghysels F->ops->matsolve = MatMatSolve_STRUMPACK; 2877d6ea485SPieter Ghysels PetscFunctionReturn(0); 2887d6ea485SPieter Ghysels } 2897d6ea485SPieter Ghysels 2907d6ea485SPieter Ghysels #undef __FUNCT__ 2917d6ea485SPieter Ghysels #define __FUNCT__ "MatFactorGetSolverPackage_aij_strumpack" 292ad0c5e61SPieter Ghysels static PetscErrorCode MatFactorGetSolverPackage_aij_strumpack(Mat A,const MatSolverPackage *type) 2937d6ea485SPieter Ghysels { 2947d6ea485SPieter Ghysels PetscFunctionBegin; 2957d6ea485SPieter Ghysels *type = MATSOLVERSTRUMPACK; 2967d6ea485SPieter Ghysels PetscFunctionReturn(0); 2977d6ea485SPieter Ghysels } 2987d6ea485SPieter Ghysels 299575704cbSPieter Ghysels /*MC 300575704cbSPieter Ghysels MATSOLVERSSTRUMPACK = "strumpack" - A solver package providing a direct sparse solver (PCLU) 301575704cbSPieter Ghysels and a preconditioner (PCILU) using low-rank compression via the external package STRUMPACK. 302575704cbSPieter Ghysels 303575704cbSPieter Ghysels Consult the STRUMPACK-sparse manual for more info. 304575704cbSPieter Ghysels 305575704cbSPieter Ghysels Use 306575704cbSPieter Ghysels ./configure --download-strumpack 307575704cbSPieter Ghysels to have PETSc installed with STRUMPACK 308575704cbSPieter Ghysels 309575704cbSPieter Ghysels Use 310575704cbSPieter Ghysels -pc_type lu -pc_factor_mat_solver_package strumpack 311575704cbSPieter Ghysels to use this as an exact (direct) solver, use 312575704cbSPieter Ghysels -pc_type ilu -pc_factor_mat_solver_package strumpack 313575704cbSPieter Ghysels to enable low-rank compression (i.e, use as a preconditioner). 314575704cbSPieter Ghysels 315575704cbSPieter Ghysels Works with AIJ matrices 316575704cbSPieter Ghysels 317575704cbSPieter Ghysels Options Database Keys: 318575704cbSPieter Ghysels + -mat_strumpack_rctol <1e-4> - Relative compression tolerance (None) 319575704cbSPieter Ghysels . -mat_strumpack_hssminsize <512> - Minimum size of dense block for HSS compression (None) 320575704cbSPieter Ghysels . -mat_strumpack_colperm <TRUE> - Permute matrix to make diagonal nonzeros (None) 321575704cbSPieter Ghysels 322575704cbSPieter Ghysels Level: beginner 323575704cbSPieter Ghysels 324575704cbSPieter Ghysels .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, PCFactorSetMatSolverPackage(), MatSolverPackage 325575704cbSPieter Ghysels M*/ 3267d6ea485SPieter Ghysels #undef __FUNCT__ 3277d6ea485SPieter Ghysels #define __FUNCT__ "MatGetFactor_aij_strumpack" 328ad0c5e61SPieter Ghysels static PetscErrorCode MatGetFactor_aij_strumpack(Mat A,MatFactorType ftype,Mat *F) 3297d6ea485SPieter Ghysels { 3307d6ea485SPieter Ghysels Mat B; 331*0d08a34dSPieter Ghysels STRUMPACK_SparseSolver *S; 3327d6ea485SPieter Ghysels PetscErrorCode ierr; 3337d6ea485SPieter Ghysels PetscInt M=A->rmap->N,N=A->cmap->N; 3347d6ea485SPieter Ghysels STRUMPACK_INTERFACE iface; 335575704cbSPieter Ghysels PetscBool verb,flg,set; 336575704cbSPieter Ghysels PetscReal rctol; 337575704cbSPieter Ghysels PetscInt hssminsize; 3387d6ea485SPieter Ghysels int argc; 33955c022e5SPieter Ghysels char **args,*copts,*pname; 34055c022e5SPieter Ghysels size_t len; 3417d6ea485SPieter Ghysels const STRUMPACK_PRECISION table[2][2][2] = {{{STRUMPACK_FLOATCOMPLEX_64,STRUMPACK_DOUBLECOMPLEX_64}, 3427d6ea485SPieter Ghysels {STRUMPACK_FLOAT_64,STRUMPACK_DOUBLE_64}}, 3437d6ea485SPieter Ghysels {{STRUMPACK_FLOATCOMPLEX,STRUMPACK_DOUBLECOMPLEX}, 3447d6ea485SPieter Ghysels {STRUMPACK_FLOAT,STRUMPACK_DOUBLE}}}; 3457d6ea485SPieter Ghysels const STRUMPACK_PRECISION prec = table[(sizeof(PetscInt)==8)?0:1][(PETSC_SCALAR==PETSC_COMPLEX)?0:1][(PETSC_REAL==PETSC_FLOAT)?0:1]; 3467d6ea485SPieter Ghysels 3477d6ea485SPieter Ghysels PetscFunctionBegin; 3487d6ea485SPieter Ghysels /* Create the factorization matrix */ 3497d6ea485SPieter Ghysels ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 3507d6ea485SPieter Ghysels ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,M,N);CHKERRQ(ierr); 3517d6ea485SPieter Ghysels ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 3527d6ea485SPieter Ghysels ierr = MatSeqAIJSetPreallocation(B,0,NULL); 3537d6ea485SPieter Ghysels ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr); 354575704cbSPieter Ghysels if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) { 3557d6ea485SPieter Ghysels B->ops->lufactorsymbolic = MatLUFactorSymbolic_STRUMPACK; 356575704cbSPieter Ghysels B->ops->ilufactorsymbolic = MatLUFactorSymbolic_STRUMPACK; 357575704cbSPieter Ghysels } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 3587d6ea485SPieter Ghysels B->ops->view = MatView_STRUMPACK; 3597d6ea485SPieter Ghysels B->ops->destroy = MatDestroy_STRUMPACK; 3607d6ea485SPieter Ghysels B->ops->getdiagonal = MatGetDiagonal_STRUMPACK; 3617d6ea485SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_aij_strumpack);CHKERRQ(ierr); 362575704cbSPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSRelCompTol_C",MatSTRUMPACKSetHSSRelCompTol_STRUMPACK);CHKERRQ(ierr); 363575704cbSPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSMinSize_C",MatSTRUMPACKSetHSSMinSize_STRUMPACK);CHKERRQ(ierr); 364575704cbSPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetColPerm_C",MatSTRUMPACKSetColPerm_STRUMPACK);CHKERRQ(ierr); 365575704cbSPieter Ghysels B->factortype = ftype; 366*0d08a34dSPieter Ghysels ierr = PetscNewLog(B,&S);CHKERRQ(ierr); 367*0d08a34dSPieter Ghysels B->spptr = S; 368*0d08a34dSPieter Ghysels 369*0d08a34dSPieter Ghysels ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg); 370*0d08a34dSPieter Ghysels iface = flg ? STRUMPACK_MT : STRUMPACK_MPI_DIST; 3717d6ea485SPieter Ghysels 3727d6ea485SPieter Ghysels ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"STRUMPACK Options","Mat");CHKERRQ(ierr); 3737d6ea485SPieter Ghysels 374575704cbSPieter Ghysels verb = PetscLogPrintInfo ? PETSC_TRUE : PETSC_FALSE; 3757d6ea485SPieter Ghysels ierr = PetscOptionsBool("-mat_strumpack_verbose","Print STRUMPACK information","None",verb,&verb,NULL);CHKERRQ(ierr); 3767d6ea485SPieter Ghysels 37755c022e5SPieter Ghysels ierr = PetscOptionsGetAll(NULL,&copts);CHKERRQ(ierr); 37855c022e5SPieter Ghysels ierr = PetscStrlen(copts,&len);CHKERRQ(ierr); 37955c022e5SPieter Ghysels len += PETSC_MAX_PATH_LEN+1; 38055c022e5SPieter Ghysels ierr = PetscMalloc1(len,&pname);CHKERRQ(ierr); 38155c022e5SPieter Ghysels /* first string is assumed to be the program name, so add program name to options string */ 38255c022e5SPieter Ghysels ierr = PetscGetProgramName(pname,len);CHKERRQ(ierr); 38355c022e5SPieter Ghysels ierr = PetscStrcat(pname," ");CHKERRQ(ierr); 38455c022e5SPieter Ghysels ierr = PetscStrcat(pname,copts);CHKERRQ(ierr); 385ee5a0f3aSPieter Ghysels ierr = PetscStrToArray(pname,' ',&argc,&args);CHKERRQ(ierr); 38655c022e5SPieter Ghysels ierr = PetscFree(copts);CHKERRQ(ierr); 38755c022e5SPieter Ghysels ierr = PetscFree(pname);CHKERRQ(ierr); 38855c022e5SPieter Ghysels 389*0d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_init",STRUMPACK_init(S,PetscObjectComm((PetscObject)A),prec,iface,argc,args,verb)); 3907d6ea485SPieter Ghysels 391*0d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_get_rctol",rctol = (PetscReal)STRUMPACK_get_rctol(*S)); 392575704cbSPieter Ghysels ierr = PetscOptionsReal("-mat_strumpack_rctol","Relative compression tolerance","None",rctol,&rctol,&set);CHKERRQ(ierr); 393*0d08a34dSPieter Ghysels if (set) PetscStackCall("STRUMPACK_set_rctol",STRUMPACK_set_rctol(*S,(double)rctol)); 394575704cbSPieter Ghysels 395*0d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_get_mc64job",flg = (STRUMPACK_get_mc64job(*S) == 0) ? PETSC_FALSE : PETSC_TRUE); 396575704cbSPieter Ghysels ierr = PetscOptionsBool("-mat_strumpack_colperm","Find a col perm to get nonzero diagonal","None",flg,&flg,&set);CHKERRQ(ierr); 397*0d08a34dSPieter Ghysels if (set) PetscStackCall("STRUMPACK_set_mc64job",STRUMPACK_set_mc64job(*S,flg ? 5 : 0)); 398575704cbSPieter Ghysels 399*0d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_get_minimum_HSS_size",hssminsize = (PetscInt)STRUMPACK_get_minimum_HSS_size(*S)); 400575704cbSPieter Ghysels ierr = PetscOptionsInt("-mat_strumpack_hssminsize","Minimum size of dense block for HSS compression","None",hssminsize,&hssminsize,&set);CHKERRQ(ierr); 401*0d08a34dSPieter Ghysels if (set) PetscStackCall("STRUMPACK_set_minimum_HSS_size",STRUMPACK_set_minimum_HSS_size(*S,(int)hssminsize)); 402575704cbSPieter Ghysels 403575704cbSPieter Ghysels PetscOptionsEnd(); 404575704cbSPieter Ghysels 40588382b05SPieter Ghysels if (ftype == MAT_FACTOR_ILU) { 40688382b05SPieter Ghysels /* When enabling HSS compression, the STRUMPACK solver becomes an incomplete */ 40788382b05SPieter Ghysels /* (or approximate) LU factorization. */ 408*0d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_use_HSS",STRUMPACK_use_HSS(*S,1)); 409575704cbSPieter Ghysels /* Disable the outer iterative solver from STRUMPACK. */ 410575704cbSPieter Ghysels /* When STRUMPACK is used as a direct solver, it will by default do iterative refinement. */ 411575704cbSPieter Ghysels /* When STRUMPACK is used with as an approximate factorization preconditioner (by enabling */ 41288382b05SPieter Ghysels /* low-rank compression), it will use it's own GMRES. Here we can disable the */ 413575704cbSPieter Ghysels /* outer iterative solver, as PETSc uses STRUMPACK from within a KSP. */ 414*0d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_set_Krylov_solver", STRUMPACK_set_Krylov_solver(*S, STRUMPACK_DIRECT)); 41588382b05SPieter Ghysels } 416575704cbSPieter Ghysels 417575704cbSPieter Ghysels /* Allow the user to set or overwrite the above options from the command line */ 418*0d08a34dSPieter Ghysels PetscStackCall("STRUMPACK_set_from_options",STRUMPACK_set_from_options(*S)); 41955c022e5SPieter Ghysels ierr = PetscStrToArrayDestroy(argc,args);CHKERRQ(ierr); 42055c022e5SPieter Ghysels 4217d6ea485SPieter Ghysels *F = B; 4227d6ea485SPieter Ghysels PetscFunctionReturn(0); 4237d6ea485SPieter Ghysels } 4247d6ea485SPieter Ghysels 4257d6ea485SPieter Ghysels #undef __FUNCT__ 4267d6ea485SPieter Ghysels #define __FUNCT__ "MatSolverPackageRegister_STRUMPACK" 4277d6ea485SPieter Ghysels PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_STRUMPACK(void) 4287d6ea485SPieter Ghysels { 4297d6ea485SPieter Ghysels PetscErrorCode ierr; 4307d6ea485SPieter Ghysels 4317d6ea485SPieter Ghysels PetscFunctionBegin; 4327d6ea485SPieter Ghysels ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr); 4337d6ea485SPieter Ghysels ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr); 434575704cbSPieter Ghysels ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_ILU,MatGetFactor_aij_strumpack);CHKERRQ(ierr); 435575704cbSPieter Ghysels ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_ILU,MatGetFactor_aij_strumpack);CHKERRQ(ierr); 4367d6ea485SPieter Ghysels PetscFunctionReturn(0); 4377d6ea485SPieter Ghysels } 438