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 /* 67d6ea485SPieter Ghysels These are only relevant for MATMPIAIJ, not for MATSEQAIJ. 77d6ea485SPieter Ghysels REPLICATED - STRUMPACK expects the entire sparse matrix and right-hand side on every process. 87d6ea485SPieter Ghysels DISTRIBUTED - STRUMPACK expects the sparse matrix and right-hand side to be distributed across the entire MPI communicator. 97d6ea485SPieter Ghysels */ 107d6ea485SPieter Ghysels typedef enum {REPLICATED, DISTRIBUTED} STRUMPACK_MatInputMode; 117d6ea485SPieter Ghysels const char *STRUMPACK_MatInputModes[] = {"REPLICATED","DISTRIBUTED","STRUMPACK_MatInputMode","PETSC_",0}; 127d6ea485SPieter Ghysels 137d6ea485SPieter Ghysels typedef struct { 147d6ea485SPieter Ghysels STRUMPACK_SparseSolver S; 157d6ea485SPieter Ghysels STRUMPACK_MatInputMode MatInputMode; 167d6ea485SPieter Ghysels } STRUMPACK_data; 177d6ea485SPieter Ghysels 187d6ea485SPieter Ghysels 197d6ea485SPieter Ghysels #undef __FUNCT__ 207d6ea485SPieter Ghysels #define __FUNCT__ "MatGetDiagonal_STRUMPACK" 21*ad0c5e61SPieter Ghysels static PetscErrorCode MatGetDiagonal_STRUMPACK(Mat A,Vec v) 227d6ea485SPieter Ghysels { 237d6ea485SPieter Ghysels PetscFunctionBegin; 247d6ea485SPieter Ghysels SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type: STRUMPACK factor"); 257d6ea485SPieter Ghysels PetscFunctionReturn(0); 267d6ea485SPieter Ghysels } 277d6ea485SPieter Ghysels 287d6ea485SPieter Ghysels #undef __FUNCT__ 297d6ea485SPieter Ghysels #define __FUNCT__ "MatDestroy_STRUMPACK" 30*ad0c5e61SPieter Ghysels static PetscErrorCode MatDestroy_STRUMPACK(Mat A) 317d6ea485SPieter Ghysels { 327d6ea485SPieter Ghysels STRUMPACK_data *sp = (STRUMPACK_data*)A->spptr; 337d6ea485SPieter Ghysels PetscErrorCode ierr; 347d6ea485SPieter Ghysels PetscBool flg; 357d6ea485SPieter Ghysels 367d6ea485SPieter Ghysels PetscFunctionBegin; 377d6ea485SPieter Ghysels /* Deallocate STRUMPACK storage */ 387d6ea485SPieter Ghysels PetscStackCall("STRUMPACK_destroy",STRUMPACK_destroy(&(sp->S))); 397d6ea485SPieter Ghysels ierr = PetscFree(A->spptr);CHKERRQ(ierr); 407d6ea485SPieter Ghysels ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr); 417d6ea485SPieter Ghysels if (flg) { 427d6ea485SPieter Ghysels ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr); 437d6ea485SPieter Ghysels } else { 447d6ea485SPieter Ghysels ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr); 457d6ea485SPieter Ghysels } 467d6ea485SPieter Ghysels PetscFunctionReturn(0); 477d6ea485SPieter Ghysels } 487d6ea485SPieter Ghysels 497d6ea485SPieter Ghysels #undef __FUNCT__ 507d6ea485SPieter Ghysels #define __FUNCT__ "MatSolve_STRUMPACK" 51*ad0c5e61SPieter Ghysels static PetscErrorCode MatSolve_STRUMPACK(Mat A,Vec b_mpi,Vec x) 527d6ea485SPieter Ghysels { 537d6ea485SPieter Ghysels STRUMPACK_data *sp = (STRUMPACK_data*)A->spptr; 547d6ea485SPieter Ghysels STRUMPACK_RETURN_CODE sp_err; 557d6ea485SPieter Ghysels PetscErrorCode ierr; 567d6ea485SPieter Ghysels PetscMPIInt size; 577d6ea485SPieter Ghysels PetscInt N=A->cmap->N; 587d6ea485SPieter Ghysels const PetscScalar *bptr; 597d6ea485SPieter Ghysels PetscScalar *xptr; 607d6ea485SPieter Ghysels Vec x_seq,b_seq; 617d6ea485SPieter Ghysels IS iden; 627d6ea485SPieter Ghysels VecScatter scat; 637d6ea485SPieter Ghysels 647d6ea485SPieter Ghysels PetscFunctionBegin; 657d6ea485SPieter Ghysels ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 667d6ea485SPieter Ghysels if (size > 1 && sp->MatInputMode == REPLICATED) { 677d6ea485SPieter Ghysels ierr = VecCreateSeq(PETSC_COMM_SELF,N,&x_seq);CHKERRQ(ierr); 687d6ea485SPieter Ghysels ierr = VecGetArray(x_seq,&xptr);CHKERRQ(ierr); 6924a4cbefSPieter Ghysels /* replicated mat input, convert b to b_seq */ 707d6ea485SPieter Ghysels ierr = VecCreateSeq(PETSC_COMM_SELF,N,&b_seq);CHKERRQ(ierr); 717d6ea485SPieter Ghysels ierr = ISCreateStride(PETSC_COMM_SELF,N,0,1,&iden);CHKERRQ(ierr); 727d6ea485SPieter Ghysels ierr = VecScatterCreate(b_mpi,iden,b_seq,iden,&scat);CHKERRQ(ierr); 737d6ea485SPieter Ghysels ierr = ISDestroy(&iden);CHKERRQ(ierr); 747d6ea485SPieter Ghysels ierr = VecScatterBegin(scat,b_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 757d6ea485SPieter Ghysels ierr = VecScatterEnd(scat,b_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 767d6ea485SPieter Ghysels ierr = VecGetArrayRead(b_seq,&bptr);CHKERRQ(ierr); 777d6ea485SPieter Ghysels } else { /* size==1 || distributed mat input */ 787d6ea485SPieter Ghysels ierr = VecGetArray(x,&xptr);CHKERRQ(ierr); 797d6ea485SPieter Ghysels ierr = VecGetArrayRead(b_mpi,&bptr);CHKERRQ(ierr); 807d6ea485SPieter Ghysels } 817d6ea485SPieter Ghysels 827d6ea485SPieter Ghysels PetscStackCall("STRUMPACK_solve",sp_err = STRUMPACK_solve(sp->S,(PetscScalar*)bptr,xptr,0)); 837d6ea485SPieter Ghysels 847d6ea485SPieter Ghysels if (sp_err != STRUMPACK_SUCCESS) { 857d6ea485SPieter Ghysels if (sp_err == STRUMPACK_MATRIX_NOT_SET) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); 867d6ea485SPieter Ghysels else if (sp_err == STRUMPACK_REORDERING_ERROR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); 877d6ea485SPieter Ghysels else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: solve failed"); 887d6ea485SPieter Ghysels } 897d6ea485SPieter Ghysels 907d6ea485SPieter Ghysels if (size > 1 && sp->MatInputMode == REPLICATED) { 917d6ea485SPieter Ghysels ierr = VecRestoreArrayRead(b_seq,&bptr);CHKERRQ(ierr); 927d6ea485SPieter Ghysels ierr = VecDestroy(&b_seq);CHKERRQ(ierr); 937d6ea485SPieter Ghysels /* convert seq x to mpi x */ 947d6ea485SPieter Ghysels ierr = VecRestoreArray(x_seq,&xptr);CHKERRQ(ierr); 957d6ea485SPieter Ghysels ierr = VecScatterBegin(scat,x_seq,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 967d6ea485SPieter Ghysels ierr = VecScatterEnd(scat,x_seq,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 977d6ea485SPieter Ghysels ierr = VecScatterDestroy(&scat);CHKERRQ(ierr); 987d6ea485SPieter Ghysels ierr = VecDestroy(&x_seq);CHKERRQ(ierr); 997d6ea485SPieter Ghysels } else { 1007d6ea485SPieter Ghysels ierr = VecRestoreArray(x,&xptr);CHKERRQ(ierr); 1017d6ea485SPieter Ghysels ierr = VecRestoreArrayRead(b_mpi,&bptr);CHKERRQ(ierr); 1027d6ea485SPieter Ghysels } 1037d6ea485SPieter Ghysels 1047d6ea485SPieter Ghysels PetscFunctionReturn(0); 1057d6ea485SPieter Ghysels } 1067d6ea485SPieter Ghysels 1077d6ea485SPieter Ghysels #undef __FUNCT__ 1087d6ea485SPieter Ghysels #define __FUNCT__ "MatMatSolve_STRUMPACK" 109*ad0c5e61SPieter Ghysels static PetscErrorCode MatMatSolve_STRUMPACK(Mat A,Mat B_mpi,Mat X) 1107d6ea485SPieter Ghysels { 1117d6ea485SPieter Ghysels PetscErrorCode ierr; 1127d6ea485SPieter Ghysels PetscBool flg; 1137d6ea485SPieter Ghysels 1147d6ea485SPieter Ghysels PetscFunctionBegin; 1157d6ea485SPieter Ghysels ierr = PetscObjectTypeCompareAny((PetscObject)B_mpi,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 1167d6ea485SPieter Ghysels if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 1177d6ea485SPieter Ghysels ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 1187d6ea485SPieter Ghysels if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); 1197d6ea485SPieter Ghysels SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_STRUMPACK() is not implemented yet"); 1207d6ea485SPieter Ghysels PetscFunctionReturn(0); 1217d6ea485SPieter Ghysels } 1227d6ea485SPieter Ghysels 123*ad0c5e61SPieter Ghysels #undef __FUNCT__ 124*ad0c5e61SPieter Ghysels #define __FUNCT__ "MatFactorInfo_STRUMPACK" 125*ad0c5e61SPieter Ghysels static PetscErrorCode MatFactorInfo_STRUMPACK(Mat A,PetscViewer viewer) 126*ad0c5e61SPieter Ghysels { 127*ad0c5e61SPieter Ghysels PetscErrorCode ierr; 128*ad0c5e61SPieter Ghysels 129*ad0c5e61SPieter Ghysels PetscFunctionBegin; 130*ad0c5e61SPieter Ghysels /* check if matrix is strumpack type */ 131*ad0c5e61SPieter Ghysels if (A->ops->solve != MatSolve_STRUMPACK) PetscFunctionReturn(0); 132*ad0c5e61SPieter Ghysels ierr = PetscViewerASCIIPrintf(viewer,"STRUMPACK sparse solver!\n");CHKERRQ(ierr); 133*ad0c5e61SPieter Ghysels PetscFunctionReturn(0); 134*ad0c5e61SPieter Ghysels } 135*ad0c5e61SPieter Ghysels 136*ad0c5e61SPieter Ghysels #undef __FUNCT__ 137*ad0c5e61SPieter Ghysels #define __FUNCT__ "MatView_STRUMPACK" 138*ad0c5e61SPieter Ghysels static PetscErrorCode MatView_STRUMPACK(Mat A,PetscViewer viewer) 139*ad0c5e61SPieter Ghysels { 140*ad0c5e61SPieter Ghysels PetscErrorCode ierr; 141*ad0c5e61SPieter Ghysels PetscBool iascii; 142*ad0c5e61SPieter Ghysels PetscViewerFormat format; 143*ad0c5e61SPieter Ghysels 144*ad0c5e61SPieter Ghysels PetscFunctionBegin; 145*ad0c5e61SPieter Ghysels ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 146*ad0c5e61SPieter Ghysels if (iascii) { 147*ad0c5e61SPieter Ghysels ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 148*ad0c5e61SPieter Ghysels if (format == PETSC_VIEWER_ASCII_INFO) { 149*ad0c5e61SPieter Ghysels ierr = MatFactorInfo_STRUMPACK(A,viewer);CHKERRQ(ierr); 150*ad0c5e61SPieter Ghysels } 151*ad0c5e61SPieter Ghysels } 152*ad0c5e61SPieter Ghysels PetscFunctionReturn(0); 153*ad0c5e61SPieter Ghysels } 1547d6ea485SPieter Ghysels 1557d6ea485SPieter Ghysels #undef __FUNCT__ 1567d6ea485SPieter Ghysels #define __FUNCT__ "MatLUFactorNumeric_STRUMPACK" 157*ad0c5e61SPieter Ghysels static PetscErrorCode MatLUFactorNumeric_STRUMPACK(Mat F,Mat A,const MatFactorInfo *info) 1587d6ea485SPieter Ghysels { 1597d6ea485SPieter Ghysels STRUMPACK_data *sp = (STRUMPACK_data*)F->spptr; 1607d6ea485SPieter Ghysels STRUMPACK_RETURN_CODE sp_err; 1617d6ea485SPieter Ghysels Mat *tseq,A_seq = NULL; 1627d6ea485SPieter Ghysels Mat_SeqAIJ *A_d,*A_o; 1637d6ea485SPieter Ghysels Mat_MPIAIJ *mat; 1647d6ea485SPieter Ghysels PetscErrorCode ierr; 1657d6ea485SPieter Ghysels PetscInt M=A->rmap->N,m=A->rmap->n,N=A->cmap->N; 1667d6ea485SPieter Ghysels PetscMPIInt size; 1677d6ea485SPieter Ghysels IS isrow; 1687d6ea485SPieter Ghysels PetscBool flg; 1697d6ea485SPieter Ghysels 1707d6ea485SPieter Ghysels PetscFunctionBegin; 1717d6ea485SPieter Ghysels ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 1727d6ea485SPieter Ghysels 1737d6ea485SPieter Ghysels ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); 1747d6ea485SPieter Ghysels if (flg) { /* A is MATMPIAIJ */ 1757d6ea485SPieter Ghysels if (sp->MatInputMode == REPLICATED) { 1767d6ea485SPieter Ghysels if (size > 1) { /* convert mpi A to seq mat A */ 1777d6ea485SPieter Ghysels ierr = ISCreateStride(PETSC_COMM_SELF,M,0,1,&isrow);CHKERRQ(ierr); 1787d6ea485SPieter Ghysels ierr = MatGetSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&tseq);CHKERRQ(ierr); 1797d6ea485SPieter Ghysels ierr = ISDestroy(&isrow);CHKERRQ(ierr); 1807d6ea485SPieter Ghysels A_seq = *tseq; 1817d6ea485SPieter Ghysels ierr = PetscFree(tseq);CHKERRQ(ierr); 1827d6ea485SPieter Ghysels A_d = (Mat_SeqAIJ*)A_seq->data; 1837d6ea485SPieter Ghysels } else { /* size == 1 */ 1847d6ea485SPieter Ghysels mat = (Mat_MPIAIJ*)A->data; 1857d6ea485SPieter Ghysels A_d = (Mat_SeqAIJ*)(mat->A)->data; 1867d6ea485SPieter Ghysels } 1877d6ea485SPieter Ghysels PetscStackCall("STRUMPACK_set_csr_matrix",STRUMPACK_set_csr_matrix(sp->S,&N,A_d->i,A_d->j,A_d->a,0)); 1887d6ea485SPieter Ghysels } else { /* sp->MatInputMode == DISTRIBUTED */ 1897d6ea485SPieter Ghysels mat = (Mat_MPIAIJ*)A->data; 1907d6ea485SPieter Ghysels A_d = (Mat_SeqAIJ*)(mat->A)->data; 1917d6ea485SPieter Ghysels A_o = (Mat_SeqAIJ*)(mat->B)->data; 1927d6ea485SPieter Ghysels PetscStackCall("STRUMPACK_set_MPIAIJ_matrix",STRUMPACK_set_MPIAIJ_matrix(sp->S,&m,A_d->i,A_d->j,A_d->a,A_o->i,A_o->j,A_o->a,mat->garray)); 1937d6ea485SPieter Ghysels } 1947d6ea485SPieter Ghysels } else { /* A is MATSEQAIJ */ 1957d6ea485SPieter Ghysels A_d = (Mat_SeqAIJ*)A->data; 1967d6ea485SPieter Ghysels PetscStackCall("STRUMPACK_set_csr_matrix",STRUMPACK_set_csr_matrix(sp->S,&N,A_d->i,A_d->j,A_d->a,0)); 1977d6ea485SPieter Ghysels } 1987d6ea485SPieter Ghysels 1997d6ea485SPieter Ghysels /* Reorder and Factor the matrix. */ 2007d6ea485SPieter Ghysels /* TODO figure out how to avoid reorder if the matrix values changed, but the pattern remains the same. */ 2017d6ea485SPieter Ghysels PetscStackCall("STRUMPACK_reorder",sp_err = STRUMPACK_reorder(sp->S)); 2027d6ea485SPieter Ghysels PetscStackCall("STRUMPACK_factor",sp_err = STRUMPACK_factor(sp->S)); 2037d6ea485SPieter Ghysels if (sp_err != STRUMPACK_SUCCESS) { 2047d6ea485SPieter Ghysels if (sp_err == STRUMPACK_MATRIX_NOT_SET) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); 2057d6ea485SPieter Ghysels else if (sp_err == STRUMPACK_REORDERING_ERROR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); 2067d6ea485SPieter Ghysels else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: factorization failed"); 2077d6ea485SPieter Ghysels } 2087d6ea485SPieter Ghysels if (flg && sp->MatInputMode == REPLICATED && size > 1) { 2097d6ea485SPieter Ghysels ierr = MatDestroy(&A_seq);CHKERRQ(ierr); 2107d6ea485SPieter Ghysels } 2117d6ea485SPieter Ghysels PetscFunctionReturn(0); 2127d6ea485SPieter Ghysels } 2137d6ea485SPieter Ghysels 2147d6ea485SPieter Ghysels #undef __FUNCT__ 2157d6ea485SPieter Ghysels #define __FUNCT__ "MatLUFactorSymbolic_STRUMPACK" 216*ad0c5e61SPieter Ghysels static PetscErrorCode MatLUFactorSymbolic_STRUMPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 2177d6ea485SPieter Ghysels { 2187d6ea485SPieter Ghysels PetscFunctionBegin; 2197d6ea485SPieter Ghysels F->ops->lufactornumeric = MatLUFactorNumeric_STRUMPACK; 2207d6ea485SPieter Ghysels F->ops->solve = MatSolve_STRUMPACK; 2217d6ea485SPieter Ghysels F->ops->matsolve = MatMatSolve_STRUMPACK; 2227d6ea485SPieter Ghysels PetscFunctionReturn(0); 2237d6ea485SPieter Ghysels } 2247d6ea485SPieter Ghysels 2257d6ea485SPieter Ghysels #undef __FUNCT__ 2267d6ea485SPieter Ghysels #define __FUNCT__ "MatFactorGetSolverPackage_aij_strumpack" 227*ad0c5e61SPieter Ghysels static PetscErrorCode MatFactorGetSolverPackage_aij_strumpack(Mat A,const MatSolverPackage *type) 2287d6ea485SPieter Ghysels { 2297d6ea485SPieter Ghysels PetscFunctionBegin; 2307d6ea485SPieter Ghysels *type = MATSOLVERSTRUMPACK; 2317d6ea485SPieter Ghysels PetscFunctionReturn(0); 2327d6ea485SPieter Ghysels } 2337d6ea485SPieter Ghysels 2347d6ea485SPieter Ghysels #undef __FUNCT__ 2357d6ea485SPieter Ghysels #define __FUNCT__ "MatGetFactor_aij_strumpack" 236*ad0c5e61SPieter Ghysels static PetscErrorCode MatGetFactor_aij_strumpack(Mat A,MatFactorType ftype,Mat *F) 2377d6ea485SPieter Ghysels { 2387d6ea485SPieter Ghysels Mat B; 2397d6ea485SPieter Ghysels STRUMPACK_data *sp; 2407d6ea485SPieter Ghysels PetscErrorCode ierr; 2417d6ea485SPieter Ghysels PetscInt M=A->rmap->N,N=A->cmap->N; 2427d6ea485SPieter Ghysels PetscMPIInt size; 2437d6ea485SPieter Ghysels STRUMPACK_INTERFACE iface; 2447d6ea485SPieter Ghysels PetscBool verb,flg; 2457d6ea485SPieter Ghysels int argc; 2467d6ea485SPieter Ghysels char **args; 2477d6ea485SPieter Ghysels const STRUMPACK_PRECISION table[2][2][2] = {{{STRUMPACK_FLOATCOMPLEX_64,STRUMPACK_DOUBLECOMPLEX_64}, 2487d6ea485SPieter Ghysels {STRUMPACK_FLOAT_64,STRUMPACK_DOUBLE_64}}, 2497d6ea485SPieter Ghysels {{STRUMPACK_FLOATCOMPLEX,STRUMPACK_DOUBLECOMPLEX}, 2507d6ea485SPieter Ghysels {STRUMPACK_FLOAT,STRUMPACK_DOUBLE}}}; 2517d6ea485SPieter Ghysels const STRUMPACK_PRECISION prec = table[(sizeof(PetscInt)==8)?0:1][(PETSC_SCALAR==PETSC_COMPLEX)?0:1][(PETSC_REAL==PETSC_FLOAT)?0:1]; 2527d6ea485SPieter Ghysels 2537d6ea485SPieter Ghysels PetscFunctionBegin; 2547d6ea485SPieter Ghysels ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 2557d6ea485SPieter Ghysels /* Create the factorization matrix */ 2567d6ea485SPieter Ghysels ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 2577d6ea485SPieter Ghysels ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,M,N);CHKERRQ(ierr); 2587d6ea485SPieter Ghysels ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 2597d6ea485SPieter Ghysels ierr = MatSeqAIJSetPreallocation(B,0,NULL); 2607d6ea485SPieter Ghysels ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr); 2617d6ea485SPieter Ghysels B->ops->lufactorsymbolic = MatLUFactorSymbolic_STRUMPACK; 2627d6ea485SPieter Ghysels B->ops->view = MatView_STRUMPACK; 2637d6ea485SPieter Ghysels B->ops->destroy = MatDestroy_STRUMPACK; 2647d6ea485SPieter Ghysels B->ops->getdiagonal = MatGetDiagonal_STRUMPACK; 2657d6ea485SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_aij_strumpack);CHKERRQ(ierr); 2667d6ea485SPieter Ghysels B->factortype = MAT_FACTOR_LU; 2677d6ea485SPieter Ghysels ierr = PetscNewLog(B,&sp);CHKERRQ(ierr); 2687d6ea485SPieter Ghysels B->spptr = sp; 2697d6ea485SPieter Ghysels 2707d6ea485SPieter Ghysels ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"STRUMPACK Options","Mat");CHKERRQ(ierr); 2717d6ea485SPieter Ghysels sp->MatInputMode = DISTRIBUTED; 2727d6ea485SPieter Ghysels iface = STRUMPACK_MPI_DIST; 27324a4cbefSPieter Ghysels ierr = PetscOptionsEnum("-mat_strumpack_matinput","Matrix input mode (replicated or distributed)","None",STRUMPACK_MatInputModes, 2747d6ea485SPieter Ghysels (PetscEnum)sp->MatInputMode,(PetscEnum*)&sp->MatInputMode,NULL);CHKERRQ(ierr); 2757d6ea485SPieter Ghysels if (sp->MatInputMode == DISTRIBUTED && size == 1) sp->MatInputMode = REPLICATED; 2767d6ea485SPieter Ghysels if (sp->MatInputMode == DISTRIBUTED) iface = STRUMPACK_MPI_DIST; 2777d6ea485SPieter Ghysels else if (sp->MatInputMode == REPLICATED) iface = STRUMPACK_MPI; 2787d6ea485SPieter Ghysels 2797d6ea485SPieter Ghysels ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg); 2807d6ea485SPieter Ghysels if (flg) iface = STRUMPACK_MT; 2817d6ea485SPieter Ghysels 2827d6ea485SPieter Ghysels if (PetscLogPrintInfo) verb = PETSC_TRUE; 2837d6ea485SPieter Ghysels else verb = PETSC_FALSE; 2847d6ea485SPieter Ghysels ierr = PetscOptionsBool("-mat_strumpack_verbose","Print STRUMPACK information","None",verb,&verb,NULL);CHKERRQ(ierr); 2857d6ea485SPieter Ghysels PetscOptionsEnd(); 2867d6ea485SPieter Ghysels 2877d6ea485SPieter Ghysels ierr = PetscGetArgs(&argc,&args);CHKERRQ(ierr); 2887d6ea485SPieter Ghysels PetscStackCall("STRUMPACK_init",STRUMPACK_init(&(sp->S),PetscObjectComm((PetscObject)A),prec,iface,argc,args,verb)); 2897d6ea485SPieter Ghysels PetscStackCall("STRUMPACK_set_from_options",STRUMPACK_set_from_options(sp->S)); 2907d6ea485SPieter Ghysels 2917d6ea485SPieter Ghysels *F = B; 2927d6ea485SPieter Ghysels PetscFunctionReturn(0); 2937d6ea485SPieter Ghysels } 2947d6ea485SPieter Ghysels 2957d6ea485SPieter Ghysels #undef __FUNCT__ 2967d6ea485SPieter Ghysels #define __FUNCT__ "MatSolverPackageRegister_STRUMPACK" 2977d6ea485SPieter Ghysels PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_STRUMPACK(void) 2987d6ea485SPieter Ghysels { 2997d6ea485SPieter Ghysels PetscErrorCode ierr; 3007d6ea485SPieter Ghysels 3017d6ea485SPieter Ghysels PetscFunctionBegin; 3027d6ea485SPieter Ghysels ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr); 3037d6ea485SPieter Ghysels ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr); 3047d6ea485SPieter Ghysels PetscFunctionReturn(0); 3057d6ea485SPieter Ghysels } 3067d6ea485SPieter Ghysels 3077d6ea485SPieter Ghysels /*MC 3087d6ea485SPieter Ghysels MATSOLVERSSTRUMPACK - Parallel direct solver package for LU factorization 3097d6ea485SPieter Ghysels 3107d6ea485SPieter Ghysels Use ./configure --download-strumpack to have PETSc installed with STRUMPACK 3117d6ea485SPieter Ghysels 3127d6ea485SPieter Ghysels Use -pc_type lu -pc_factor_mat_solver_package strumpack to us this direct solver 3137d6ea485SPieter Ghysels 3147d6ea485SPieter Ghysels Works with AIJ matrices 3157d6ea485SPieter Ghysels 3167d6ea485SPieter Ghysels .seealso: PCLU 3177d6ea485SPieter Ghysels 3187d6ea485SPieter Ghysels .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage 3197d6ea485SPieter Ghysels 3207d6ea485SPieter Ghysels M*/ 3217d6ea485SPieter Ghysels 322