1*7d6ea485SPieter Ghysels #include <../src/mat/impls/aij/seq/aij.h> 2*7d6ea485SPieter Ghysels #include <../src/mat/impls/aij/mpi/mpiaij.h> 3*7d6ea485SPieter Ghysels #include <StrumpackSparseSolver.h> 4*7d6ea485SPieter Ghysels 5*7d6ea485SPieter Ghysels /* 6*7d6ea485SPieter Ghysels These are only relevant for MATMPIAIJ, not for MATSEQAIJ. 7*7d6ea485SPieter Ghysels REPLICATED - STRUMPACK expects the entire sparse matrix and right-hand side on every process. 8*7d6ea485SPieter Ghysels DISTRIBUTED - STRUMPACK expects the sparse matrix and right-hand side to be distributed across the entire MPI communicator. 9*7d6ea485SPieter Ghysels */ 10*7d6ea485SPieter Ghysels typedef enum {REPLICATED, DISTRIBUTED} STRUMPACK_MatInputMode; 11*7d6ea485SPieter Ghysels const char *STRUMPACK_MatInputModes[] = {"REPLICATED","DISTRIBUTED","STRUMPACK_MatInputMode","PETSC_",0}; 12*7d6ea485SPieter Ghysels 13*7d6ea485SPieter Ghysels typedef struct { 14*7d6ea485SPieter Ghysels STRUMPACK_SparseSolver S; 15*7d6ea485SPieter Ghysels STRUMPACK_MatInputMode MatInputMode; 16*7d6ea485SPieter Ghysels } STRUMPACK_data; 17*7d6ea485SPieter Ghysels 18*7d6ea485SPieter Ghysels extern PetscErrorCode MatFactorInfo_STRUMPACK(Mat,PetscViewer); 19*7d6ea485SPieter Ghysels extern PetscErrorCode MatLUFactorNumeric_STRUMPACK(Mat,Mat,const MatFactorInfo*); 20*7d6ea485SPieter Ghysels extern PetscErrorCode MatDestroy_STRUMPACK(Mat); 21*7d6ea485SPieter Ghysels extern PetscErrorCode MatView_STRUMPACK(Mat,PetscViewer); 22*7d6ea485SPieter Ghysels extern PetscErrorCode MatSolve_STRUMPACK(Mat,Vec,Vec); 23*7d6ea485SPieter Ghysels extern PetscErrorCode MatLUFactorSymbolic_STRUMPACK(Mat,Mat,IS,IS,const MatFactorInfo*); 24*7d6ea485SPieter Ghysels extern PetscErrorCode MatDestroy_MPIAIJ(Mat); 25*7d6ea485SPieter Ghysels 26*7d6ea485SPieter Ghysels #undef __FUNCT__ 27*7d6ea485SPieter Ghysels #define __FUNCT__ "MatGetDiagonal_STRUMPACK" 28*7d6ea485SPieter Ghysels PetscErrorCode MatGetDiagonal_STRUMPACK(Mat A,Vec v) 29*7d6ea485SPieter Ghysels { 30*7d6ea485SPieter Ghysels PetscFunctionBegin; 31*7d6ea485SPieter Ghysels SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type: STRUMPACK factor"); 32*7d6ea485SPieter Ghysels PetscFunctionReturn(0); 33*7d6ea485SPieter Ghysels } 34*7d6ea485SPieter Ghysels 35*7d6ea485SPieter Ghysels #undef __FUNCT__ 36*7d6ea485SPieter Ghysels #define __FUNCT__ "MatDestroy_STRUMPACK" 37*7d6ea485SPieter Ghysels PetscErrorCode MatDestroy_STRUMPACK(Mat A) 38*7d6ea485SPieter Ghysels { 39*7d6ea485SPieter Ghysels STRUMPACK_data *sp = (STRUMPACK_data*)A->spptr; 40*7d6ea485SPieter Ghysels PetscErrorCode ierr; 41*7d6ea485SPieter Ghysels PetscBool flg; 42*7d6ea485SPieter Ghysels 43*7d6ea485SPieter Ghysels PetscFunctionBegin; 44*7d6ea485SPieter Ghysels /* Deallocate STRUMPACK storage */ 45*7d6ea485SPieter Ghysels PetscStackCall("STRUMPACK_destroy",STRUMPACK_destroy(&(sp->S))); 46*7d6ea485SPieter Ghysels ierr = PetscFree(A->spptr);CHKERRQ(ierr); 47*7d6ea485SPieter Ghysels ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr); 48*7d6ea485SPieter Ghysels if (flg) { 49*7d6ea485SPieter Ghysels ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr); 50*7d6ea485SPieter Ghysels } else { 51*7d6ea485SPieter Ghysels ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr); 52*7d6ea485SPieter Ghysels } 53*7d6ea485SPieter Ghysels PetscFunctionReturn(0); 54*7d6ea485SPieter Ghysels } 55*7d6ea485SPieter Ghysels 56*7d6ea485SPieter Ghysels #undef __FUNCT__ 57*7d6ea485SPieter Ghysels #define __FUNCT__ "MatSolve_STRUMPACK" 58*7d6ea485SPieter Ghysels PetscErrorCode MatSolve_STRUMPACK(Mat A,Vec b_mpi,Vec x) 59*7d6ea485SPieter Ghysels { 60*7d6ea485SPieter Ghysels STRUMPACK_data *sp = (STRUMPACK_data*)A->spptr; 61*7d6ea485SPieter Ghysels STRUMPACK_RETURN_CODE sp_err; 62*7d6ea485SPieter Ghysels PetscErrorCode ierr; 63*7d6ea485SPieter Ghysels PetscMPIInt size; 64*7d6ea485SPieter Ghysels PetscInt N=A->cmap->N; 65*7d6ea485SPieter Ghysels const PetscScalar *bptr; 66*7d6ea485SPieter Ghysels PetscScalar *xptr; 67*7d6ea485SPieter Ghysels Vec x_seq,b_seq; 68*7d6ea485SPieter Ghysels IS iden; 69*7d6ea485SPieter Ghysels VecScatter scat; 70*7d6ea485SPieter Ghysels 71*7d6ea485SPieter Ghysels PetscFunctionBegin; 72*7d6ea485SPieter Ghysels ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 73*7d6ea485SPieter Ghysels if (size > 1 && sp->MatInputMode == REPLICATED) { 74*7d6ea485SPieter Ghysels ierr = VecCreateSeq(PETSC_COMM_SELF,N,&x_seq);CHKERRQ(ierr); 75*7d6ea485SPieter Ghysels ierr = VecGetArray(x_seq,&xptr);CHKERRQ(ierr); 76*7d6ea485SPieter Ghysels /* global mat input, convert b to b_seq */ 77*7d6ea485SPieter Ghysels ierr = VecCreateSeq(PETSC_COMM_SELF,N,&b_seq);CHKERRQ(ierr); 78*7d6ea485SPieter Ghysels ierr = ISCreateStride(PETSC_COMM_SELF,N,0,1,&iden);CHKERRQ(ierr); 79*7d6ea485SPieter Ghysels ierr = VecScatterCreate(b_mpi,iden,b_seq,iden,&scat);CHKERRQ(ierr); 80*7d6ea485SPieter Ghysels ierr = ISDestroy(&iden);CHKERRQ(ierr); 81*7d6ea485SPieter Ghysels ierr = VecScatterBegin(scat,b_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 82*7d6ea485SPieter Ghysels ierr = VecScatterEnd(scat,b_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 83*7d6ea485SPieter Ghysels ierr = VecGetArrayRead(b_seq,&bptr);CHKERRQ(ierr); 84*7d6ea485SPieter Ghysels } else { /* size==1 || distributed mat input */ 85*7d6ea485SPieter Ghysels ierr = VecGetArray(x,&xptr);CHKERRQ(ierr); 86*7d6ea485SPieter Ghysels ierr = VecGetArrayRead(b_mpi,&bptr);CHKERRQ(ierr); 87*7d6ea485SPieter Ghysels } 88*7d6ea485SPieter Ghysels 89*7d6ea485SPieter Ghysels PetscStackCall("STRUMPACK_solve",sp_err = STRUMPACK_solve(sp->S,(PetscScalar*)bptr,xptr,0)); 90*7d6ea485SPieter Ghysels 91*7d6ea485SPieter Ghysels if (sp_err != STRUMPACK_SUCCESS) { 92*7d6ea485SPieter Ghysels if (sp_err == STRUMPACK_MATRIX_NOT_SET) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); 93*7d6ea485SPieter Ghysels else if (sp_err == STRUMPACK_REORDERING_ERROR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); 94*7d6ea485SPieter Ghysels else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: solve failed"); 95*7d6ea485SPieter Ghysels } 96*7d6ea485SPieter Ghysels 97*7d6ea485SPieter Ghysels if (size > 1 && sp->MatInputMode == REPLICATED) { 98*7d6ea485SPieter Ghysels ierr = VecRestoreArrayRead(b_seq,&bptr);CHKERRQ(ierr); 99*7d6ea485SPieter Ghysels ierr = VecDestroy(&b_seq);CHKERRQ(ierr); 100*7d6ea485SPieter Ghysels /* convert seq x to mpi x */ 101*7d6ea485SPieter Ghysels ierr = VecRestoreArray(x_seq,&xptr);CHKERRQ(ierr); 102*7d6ea485SPieter Ghysels ierr = VecScatterBegin(scat,x_seq,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 103*7d6ea485SPieter Ghysels ierr = VecScatterEnd(scat,x_seq,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 104*7d6ea485SPieter Ghysels ierr = VecScatterDestroy(&scat);CHKERRQ(ierr); 105*7d6ea485SPieter Ghysels ierr = VecDestroy(&x_seq);CHKERRQ(ierr); 106*7d6ea485SPieter Ghysels } else { 107*7d6ea485SPieter Ghysels ierr = VecRestoreArray(x,&xptr);CHKERRQ(ierr); 108*7d6ea485SPieter Ghysels ierr = VecRestoreArrayRead(b_mpi,&bptr);CHKERRQ(ierr); 109*7d6ea485SPieter Ghysels } 110*7d6ea485SPieter Ghysels 111*7d6ea485SPieter Ghysels PetscFunctionReturn(0); 112*7d6ea485SPieter Ghysels } 113*7d6ea485SPieter Ghysels 114*7d6ea485SPieter Ghysels #undef __FUNCT__ 115*7d6ea485SPieter Ghysels #define __FUNCT__ "MatMatSolve_STRUMPACK" 116*7d6ea485SPieter Ghysels PetscErrorCode MatMatSolve_STRUMPACK(Mat A,Mat B_mpi,Mat X) 117*7d6ea485SPieter Ghysels { 118*7d6ea485SPieter Ghysels PetscErrorCode ierr; 119*7d6ea485SPieter Ghysels PetscBool flg; 120*7d6ea485SPieter Ghysels 121*7d6ea485SPieter Ghysels PetscFunctionBegin; 122*7d6ea485SPieter Ghysels ierr = PetscObjectTypeCompareAny((PetscObject)B_mpi,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 123*7d6ea485SPieter Ghysels if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 124*7d6ea485SPieter Ghysels ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 125*7d6ea485SPieter Ghysels if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); 126*7d6ea485SPieter Ghysels SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_STRUMPACK() is not implemented yet"); 127*7d6ea485SPieter Ghysels PetscFunctionReturn(0); 128*7d6ea485SPieter Ghysels } 129*7d6ea485SPieter Ghysels 130*7d6ea485SPieter Ghysels 131*7d6ea485SPieter Ghysels #undef __FUNCT__ 132*7d6ea485SPieter Ghysels #define __FUNCT__ "MatLUFactorNumeric_STRUMPACK" 133*7d6ea485SPieter Ghysels PetscErrorCode MatLUFactorNumeric_STRUMPACK(Mat F,Mat A,const MatFactorInfo *info) 134*7d6ea485SPieter Ghysels { 135*7d6ea485SPieter Ghysels STRUMPACK_data *sp = (STRUMPACK_data*)F->spptr; 136*7d6ea485SPieter Ghysels STRUMPACK_RETURN_CODE sp_err; 137*7d6ea485SPieter Ghysels Mat *tseq,A_seq = NULL; 138*7d6ea485SPieter Ghysels Mat_SeqAIJ *A_d,*A_o; 139*7d6ea485SPieter Ghysels Mat_MPIAIJ *mat; 140*7d6ea485SPieter Ghysels PetscErrorCode ierr; 141*7d6ea485SPieter Ghysels PetscInt M=A->rmap->N,m=A->rmap->n,N=A->cmap->N; 142*7d6ea485SPieter Ghysels PetscMPIInt size; 143*7d6ea485SPieter Ghysels IS isrow; 144*7d6ea485SPieter Ghysels PetscBool flg; 145*7d6ea485SPieter Ghysels 146*7d6ea485SPieter Ghysels PetscFunctionBegin; 147*7d6ea485SPieter Ghysels ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 148*7d6ea485SPieter Ghysels 149*7d6ea485SPieter Ghysels ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); 150*7d6ea485SPieter Ghysels if (flg) { /* A is MATMPIAIJ */ 151*7d6ea485SPieter Ghysels if (sp->MatInputMode == REPLICATED) { 152*7d6ea485SPieter Ghysels if (size > 1) { /* convert mpi A to seq mat A */ 153*7d6ea485SPieter Ghysels ierr = ISCreateStride(PETSC_COMM_SELF,M,0,1,&isrow);CHKERRQ(ierr); 154*7d6ea485SPieter Ghysels ierr = MatGetSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&tseq);CHKERRQ(ierr); 155*7d6ea485SPieter Ghysels ierr = ISDestroy(&isrow);CHKERRQ(ierr); 156*7d6ea485SPieter Ghysels A_seq = *tseq; 157*7d6ea485SPieter Ghysels ierr = PetscFree(tseq);CHKERRQ(ierr); 158*7d6ea485SPieter Ghysels A_d = (Mat_SeqAIJ*)A_seq->data; 159*7d6ea485SPieter Ghysels } else { /* size == 1 */ 160*7d6ea485SPieter Ghysels mat = (Mat_MPIAIJ*)A->data; 161*7d6ea485SPieter Ghysels A_d = (Mat_SeqAIJ*)(mat->A)->data; 162*7d6ea485SPieter Ghysels } 163*7d6ea485SPieter Ghysels PetscStackCall("STRUMPACK_set_csr_matrix",STRUMPACK_set_csr_matrix(sp->S,&N,A_d->i,A_d->j,A_d->a,0)); 164*7d6ea485SPieter Ghysels } else { /* sp->MatInputMode == DISTRIBUTED */ 165*7d6ea485SPieter Ghysels mat = (Mat_MPIAIJ*)A->data; 166*7d6ea485SPieter Ghysels A_d = (Mat_SeqAIJ*)(mat->A)->data; 167*7d6ea485SPieter Ghysels A_o = (Mat_SeqAIJ*)(mat->B)->data; 168*7d6ea485SPieter 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)); 169*7d6ea485SPieter Ghysels } 170*7d6ea485SPieter Ghysels } else { /* A is MATSEQAIJ */ 171*7d6ea485SPieter Ghysels A_d = (Mat_SeqAIJ*)A->data; 172*7d6ea485SPieter Ghysels PetscStackCall("STRUMPACK_set_csr_matrix",STRUMPACK_set_csr_matrix(sp->S,&N,A_d->i,A_d->j,A_d->a,0)); 173*7d6ea485SPieter Ghysels } 174*7d6ea485SPieter Ghysels 175*7d6ea485SPieter Ghysels /* Reorder and Factor the matrix. */ 176*7d6ea485SPieter Ghysels /* TODO figure out how to avoid reorder if the matrix values changed, but the pattern remains the same. */ 177*7d6ea485SPieter Ghysels PetscStackCall("STRUMPACK_reorder",sp_err = STRUMPACK_reorder(sp->S)); 178*7d6ea485SPieter Ghysels PetscStackCall("STRUMPACK_factor",sp_err = STRUMPACK_factor(sp->S)); 179*7d6ea485SPieter Ghysels if (sp_err != STRUMPACK_SUCCESS) { 180*7d6ea485SPieter Ghysels if (sp_err == STRUMPACK_MATRIX_NOT_SET) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); 181*7d6ea485SPieter Ghysels else if (sp_err == STRUMPACK_REORDERING_ERROR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); 182*7d6ea485SPieter Ghysels else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: factorization failed"); 183*7d6ea485SPieter Ghysels } 184*7d6ea485SPieter Ghysels if (flg && sp->MatInputMode == REPLICATED && size > 1) { 185*7d6ea485SPieter Ghysels ierr = MatDestroy(&A_seq);CHKERRQ(ierr); 186*7d6ea485SPieter Ghysels } 187*7d6ea485SPieter Ghysels PetscFunctionReturn(0); 188*7d6ea485SPieter Ghysels } 189*7d6ea485SPieter Ghysels 190*7d6ea485SPieter Ghysels #undef __FUNCT__ 191*7d6ea485SPieter Ghysels #define __FUNCT__ "MatLUFactorSymbolic_STRUMPACK" 192*7d6ea485SPieter Ghysels PetscErrorCode MatLUFactorSymbolic_STRUMPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 193*7d6ea485SPieter Ghysels { 194*7d6ea485SPieter Ghysels PetscFunctionBegin; 195*7d6ea485SPieter Ghysels F->ops->lufactornumeric = MatLUFactorNumeric_STRUMPACK; 196*7d6ea485SPieter Ghysels F->ops->solve = MatSolve_STRUMPACK; 197*7d6ea485SPieter Ghysels F->ops->matsolve = MatMatSolve_STRUMPACK; 198*7d6ea485SPieter Ghysels PetscFunctionReturn(0); 199*7d6ea485SPieter Ghysels } 200*7d6ea485SPieter Ghysels 201*7d6ea485SPieter Ghysels #undef __FUNCT__ 202*7d6ea485SPieter Ghysels #define __FUNCT__ "MatFactorGetSolverPackage_aij_strumpack" 203*7d6ea485SPieter Ghysels PetscErrorCode MatFactorGetSolverPackage_aij_strumpack(Mat A,const MatSolverPackage *type) 204*7d6ea485SPieter Ghysels { 205*7d6ea485SPieter Ghysels PetscFunctionBegin; 206*7d6ea485SPieter Ghysels *type = MATSOLVERSTRUMPACK; 207*7d6ea485SPieter Ghysels PetscFunctionReturn(0); 208*7d6ea485SPieter Ghysels } 209*7d6ea485SPieter Ghysels 210*7d6ea485SPieter Ghysels #undef __FUNCT__ 211*7d6ea485SPieter Ghysels #define __FUNCT__ "MatGetFactor_aij_strumpack" 212*7d6ea485SPieter Ghysels PETSC_EXTERN PetscErrorCode MatGetFactor_aij_strumpack(Mat A,MatFactorType ftype,Mat *F) 213*7d6ea485SPieter Ghysels { 214*7d6ea485SPieter Ghysels Mat B; 215*7d6ea485SPieter Ghysels STRUMPACK_data *sp; 216*7d6ea485SPieter Ghysels PetscErrorCode ierr; 217*7d6ea485SPieter Ghysels PetscInt M=A->rmap->N,N=A->cmap->N; 218*7d6ea485SPieter Ghysels PetscMPIInt size; 219*7d6ea485SPieter Ghysels STRUMPACK_INTERFACE iface; 220*7d6ea485SPieter Ghysels PetscBool verb,flg; 221*7d6ea485SPieter Ghysels int argc; 222*7d6ea485SPieter Ghysels char **args; 223*7d6ea485SPieter Ghysels const STRUMPACK_PRECISION table[2][2][2] = {{{STRUMPACK_FLOATCOMPLEX_64,STRUMPACK_DOUBLECOMPLEX_64}, 224*7d6ea485SPieter Ghysels {STRUMPACK_FLOAT_64,STRUMPACK_DOUBLE_64}}, 225*7d6ea485SPieter Ghysels {{STRUMPACK_FLOATCOMPLEX,STRUMPACK_DOUBLECOMPLEX}, 226*7d6ea485SPieter Ghysels {STRUMPACK_FLOAT,STRUMPACK_DOUBLE}}}; 227*7d6ea485SPieter Ghysels const STRUMPACK_PRECISION prec = table[(sizeof(PetscInt)==8)?0:1][(PETSC_SCALAR==PETSC_COMPLEX)?0:1][(PETSC_REAL==PETSC_FLOAT)?0:1]; 228*7d6ea485SPieter Ghysels 229*7d6ea485SPieter Ghysels PetscFunctionBegin; 230*7d6ea485SPieter Ghysels ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 231*7d6ea485SPieter Ghysels /* Create the factorization matrix */ 232*7d6ea485SPieter Ghysels ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 233*7d6ea485SPieter Ghysels ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,M,N);CHKERRQ(ierr); 234*7d6ea485SPieter Ghysels ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 235*7d6ea485SPieter Ghysels ierr = MatSeqAIJSetPreallocation(B,0,NULL); 236*7d6ea485SPieter Ghysels ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr); 237*7d6ea485SPieter Ghysels B->ops->lufactorsymbolic = MatLUFactorSymbolic_STRUMPACK; 238*7d6ea485SPieter Ghysels B->ops->view = MatView_STRUMPACK; 239*7d6ea485SPieter Ghysels B->ops->destroy = MatDestroy_STRUMPACK; 240*7d6ea485SPieter Ghysels B->ops->getdiagonal = MatGetDiagonal_STRUMPACK; 241*7d6ea485SPieter Ghysels ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_aij_strumpack);CHKERRQ(ierr); 242*7d6ea485SPieter Ghysels B->factortype = MAT_FACTOR_LU; 243*7d6ea485SPieter Ghysels ierr = PetscNewLog(B,&sp);CHKERRQ(ierr); 244*7d6ea485SPieter Ghysels B->spptr = sp; 245*7d6ea485SPieter Ghysels 246*7d6ea485SPieter Ghysels ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"STRUMPACK Options","Mat");CHKERRQ(ierr); 247*7d6ea485SPieter Ghysels sp->MatInputMode = DISTRIBUTED; 248*7d6ea485SPieter Ghysels iface = STRUMPACK_MPI_DIST; 249*7d6ea485SPieter Ghysels ierr = PetscOptionsEnum("-mat_strumpack_matinput","Matrix input mode (global or distributed)","None",STRUMPACK_MatInputModes, 250*7d6ea485SPieter Ghysels (PetscEnum)sp->MatInputMode,(PetscEnum*)&sp->MatInputMode,NULL);CHKERRQ(ierr); 251*7d6ea485SPieter Ghysels if (sp->MatInputMode == DISTRIBUTED && size == 1) sp->MatInputMode = REPLICATED; 252*7d6ea485SPieter Ghysels if (sp->MatInputMode == DISTRIBUTED) iface = STRUMPACK_MPI_DIST; 253*7d6ea485SPieter Ghysels else if (sp->MatInputMode == REPLICATED) iface = STRUMPACK_MPI; 254*7d6ea485SPieter Ghysels 255*7d6ea485SPieter Ghysels ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg); 256*7d6ea485SPieter Ghysels if (flg) iface = STRUMPACK_MT; 257*7d6ea485SPieter Ghysels 258*7d6ea485SPieter Ghysels if (PetscLogPrintInfo) verb = PETSC_TRUE; 259*7d6ea485SPieter Ghysels else verb = PETSC_FALSE; 260*7d6ea485SPieter Ghysels ierr = PetscOptionsBool("-mat_strumpack_verbose","Print STRUMPACK information","None",verb,&verb,NULL);CHKERRQ(ierr); 261*7d6ea485SPieter Ghysels PetscOptionsEnd(); 262*7d6ea485SPieter Ghysels 263*7d6ea485SPieter Ghysels ierr = PetscGetArgs(&argc,&args);CHKERRQ(ierr); 264*7d6ea485SPieter Ghysels PetscStackCall("STRUMPACK_init",STRUMPACK_init(&(sp->S),PetscObjectComm((PetscObject)A),prec,iface,argc,args,verb)); 265*7d6ea485SPieter Ghysels PetscStackCall("STRUMPACK_set_from_options",STRUMPACK_set_from_options(sp->S)); 266*7d6ea485SPieter Ghysels 267*7d6ea485SPieter Ghysels *F = B; 268*7d6ea485SPieter Ghysels PetscFunctionReturn(0); 269*7d6ea485SPieter Ghysels } 270*7d6ea485SPieter Ghysels 271*7d6ea485SPieter Ghysels #undef __FUNCT__ 272*7d6ea485SPieter Ghysels #define __FUNCT__ "MatSolverPackageRegister_STRUMPACK" 273*7d6ea485SPieter Ghysels PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_STRUMPACK(void) 274*7d6ea485SPieter Ghysels { 275*7d6ea485SPieter Ghysels PetscErrorCode ierr; 276*7d6ea485SPieter Ghysels 277*7d6ea485SPieter Ghysels PetscFunctionBegin; 278*7d6ea485SPieter Ghysels ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr); 279*7d6ea485SPieter Ghysels ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr); 280*7d6ea485SPieter Ghysels PetscFunctionReturn(0); 281*7d6ea485SPieter Ghysels } 282*7d6ea485SPieter Ghysels 283*7d6ea485SPieter Ghysels #undef __FUNCT__ 284*7d6ea485SPieter Ghysels #define __FUNCT__ "MatFactorInfo_STRUMPACK" 285*7d6ea485SPieter Ghysels PetscErrorCode MatFactorInfo_STRUMPACK(Mat A,PetscViewer viewer) 286*7d6ea485SPieter Ghysels { 287*7d6ea485SPieter Ghysels PetscErrorCode ierr; 288*7d6ea485SPieter Ghysels 289*7d6ea485SPieter Ghysels PetscFunctionBegin; 290*7d6ea485SPieter Ghysels /* check if matrix is strumpack type */ 291*7d6ea485SPieter Ghysels if (A->ops->solve != MatSolve_STRUMPACK) PetscFunctionReturn(0); 292*7d6ea485SPieter Ghysels ierr = PetscViewerASCIIPrintf(viewer,"STRUMPACK sparse solver!\n");CHKERRQ(ierr); 293*7d6ea485SPieter Ghysels PetscFunctionReturn(0); 294*7d6ea485SPieter Ghysels } 295*7d6ea485SPieter Ghysels 296*7d6ea485SPieter Ghysels #undef __FUNCT__ 297*7d6ea485SPieter Ghysels #define __FUNCT__ "MatView_STRUMPACK" 298*7d6ea485SPieter Ghysels PetscErrorCode MatView_STRUMPACK(Mat A,PetscViewer viewer) 299*7d6ea485SPieter Ghysels { 300*7d6ea485SPieter Ghysels PetscErrorCode ierr; 301*7d6ea485SPieter Ghysels PetscBool iascii; 302*7d6ea485SPieter Ghysels PetscViewerFormat format; 303*7d6ea485SPieter Ghysels 304*7d6ea485SPieter Ghysels PetscFunctionBegin; 305*7d6ea485SPieter Ghysels ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 306*7d6ea485SPieter Ghysels if (iascii) { 307*7d6ea485SPieter Ghysels ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 308*7d6ea485SPieter Ghysels if (format == PETSC_VIEWER_ASCII_INFO) { 309*7d6ea485SPieter Ghysels ierr = MatFactorInfo_STRUMPACK(A,viewer);CHKERRQ(ierr); 310*7d6ea485SPieter Ghysels } 311*7d6ea485SPieter Ghysels } 312*7d6ea485SPieter Ghysels PetscFunctionReturn(0); 313*7d6ea485SPieter Ghysels } 314*7d6ea485SPieter Ghysels 315*7d6ea485SPieter Ghysels 316*7d6ea485SPieter Ghysels /*MC 317*7d6ea485SPieter Ghysels MATSOLVERSSTRUMPACK - Parallel direct solver package for LU factorization 318*7d6ea485SPieter Ghysels 319*7d6ea485SPieter Ghysels Use ./configure --download-strumpack to have PETSc installed with STRUMPACK 320*7d6ea485SPieter Ghysels 321*7d6ea485SPieter Ghysels Use -pc_type lu -pc_factor_mat_solver_package strumpack to us this direct solver 322*7d6ea485SPieter Ghysels 323*7d6ea485SPieter Ghysels Works with AIJ matrices 324*7d6ea485SPieter Ghysels 325*7d6ea485SPieter Ghysels .seealso: PCLU 326*7d6ea485SPieter Ghysels 327*7d6ea485SPieter Ghysels .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage 328*7d6ea485SPieter Ghysels 329*7d6ea485SPieter Ghysels M*/ 330*7d6ea485SPieter Ghysels 331