13bf14a46SMatthew Knepley /* 23bf14a46SMatthew Knepley Provides an interface to the PaStiX sparse solver 33bf14a46SMatthew Knepley */ 4c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> 5c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> 6c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h> 7c6db04a5SJed Brown #include <../src/mat/impls/sbaij/mpi/mpisbaij.h> 83bf14a46SMatthew Knepley 9dbbbd53dSSatish Balay #if defined(PETSC_USE_COMPLEX) 105ec454cfSSatish Balay #define _H_COMPLEX 115ec454cfSSatish Balay #endif 125ec454cfSSatish Balay 13c6db04a5SJed Brown #include <pastix.h> 143bf14a46SMatthew Knepley 15519f805aSKarl Rupp #if defined(PETSC_USE_COMPLEX) 16519f805aSKarl Rupp #if defined(PETSC_USE_REAL_SINGLE) 17688c8ee7SFlorent Pruvost #define SPM_FLTTYPE SpmComplex32 185ec454cfSSatish Balay #else 19688c8ee7SFlorent Pruvost #define SPM_FLTTYPE SpmComplex64 205ec454cfSSatish Balay #endif 21d41469e0Sxavier lacoste #else /* PETSC_USE_COMPLEX */ 22d41469e0Sxavier lacoste 23519f805aSKarl Rupp #if defined(PETSC_USE_REAL_SINGLE) 24688c8ee7SFlorent Pruvost #define SPM_FLTTYPE SpmFloat 255ec454cfSSatish Balay #else 26688c8ee7SFlorent Pruvost #define SPM_FLTTYPE SpmDouble 275ec454cfSSatish Balay #endif 28d41469e0Sxavier lacoste 29d41469e0Sxavier lacoste #endif /* PETSC_USE_COMPLEX */ 30d41469e0Sxavier lacoste 31dbbbd53dSSatish Balay typedef PetscScalar PastixScalar; 32dbbbd53dSSatish Balay 333bf14a46SMatthew Knepley typedef struct Mat_Pastix_ { 343bf14a46SMatthew Knepley pastix_data_t *pastix_data; /* Pastix data storage structure */ 35688c8ee7SFlorent Pruvost MPI_Comm comm; /* MPI Communicator used to initialize pastix */ 36688c8ee7SFlorent Pruvost spmatrix_t *spm; /* SPM matrix structure */ 37688c8ee7SFlorent Pruvost MatStructure matstruc; /* DIFFERENT_NONZERO_PATTERN if uninitilized, SAME otherwise */ 38688c8ee7SFlorent Pruvost PetscScalar *rhs; /* Right-hand-side member */ 39688c8ee7SFlorent Pruvost PetscInt rhsnbr; /* Right-hand-side number */ 40688c8ee7SFlorent Pruvost pastix_int_t iparm[IPARM_SIZE]; /* Integer parameters */ 41baed9decSBarry Smith double dparm[DPARM_SIZE]; /* Floating point parameters */ 423bf14a46SMatthew Knepley } Mat_Pastix; 433bf14a46SMatthew Knepley 443bf14a46SMatthew Knepley /* 45688c8ee7SFlorent Pruvost convert PETSc matrix to SPM structure 463bf14a46SMatthew Knepley 473bf14a46SMatthew Knepley input: 48688c8ee7SFlorent Pruvost A - matrix in aij, baij or sbaij format 49688c8ee7SFlorent Pruvost reuse - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple 50688c8ee7SFlorent Pruvost MAT_REUSE_MATRIX: only the values in v array are updated 513bf14a46SMatthew Knepley output: 52688c8ee7SFlorent Pruvost spm - The SPM built from A 533bf14a46SMatthew Knepley */ 54688c8ee7SFlorent Pruvost static PetscErrorCode MatConvertToSPM(Mat A, MatReuse reuse, Mat_Pastix *pastix) 55d71ae5a4SJacob Faibussowitsch { 56688c8ee7SFlorent Pruvost Mat A_loc, A_aij; 57688c8ee7SFlorent Pruvost const PetscInt *row, *col; 58688c8ee7SFlorent Pruvost PetscInt n, i; 59688c8ee7SFlorent Pruvost const PetscScalar *val; 60688c8ee7SFlorent Pruvost PetscBool ismpiaij, isseqaij, ismpisbaij, isseqsbaij; 61688c8ee7SFlorent Pruvost PetscBool flag; 62688c8ee7SFlorent Pruvost spmatrix_t spm2, *spm = NULL; 63688c8ee7SFlorent Pruvost int spm_err; 643bf14a46SMatthew Knepley 653bf14a46SMatthew Knepley PetscFunctionBegin; 66688c8ee7SFlorent Pruvost /* Get A datas */ 67688c8ee7SFlorent Pruvost PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATSEQAIJ, &isseqaij)); 68688c8ee7SFlorent Pruvost PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATMPIAIJ, &ismpiaij)); 69688c8ee7SFlorent Pruvost /* TODO: Block Aij should be handled with dof in spm */ 70688c8ee7SFlorent Pruvost PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATSEQSBAIJ, &isseqsbaij)); 71688c8ee7SFlorent Pruvost PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATMPISBAIJ, &ismpisbaij)); 723bf14a46SMatthew Knepley 73688c8ee7SFlorent Pruvost if (isseqsbaij || ismpisbaij) PetscCall(MatConvert(A, MATAIJ, reuse, &A_aij)); 74688c8ee7SFlorent Pruvost else A_aij = A; 75745c78f7SBarry Smith 76688c8ee7SFlorent Pruvost if (ismpiaij || ismpisbaij) PetscCall(MatMPIAIJGetLocalMat(A_aij, MAT_INITIAL_MATRIX, &A_loc)); 77688c8ee7SFlorent Pruvost else if (isseqaij || isseqsbaij) A_loc = A_aij; 78688c8ee7SFlorent Pruvost else SETERRQ(PetscObjectComm((PetscObject)A_aij), PETSC_ERR_SUP, "Not for type %s", ((PetscObject)A)->type_name); 79688c8ee7SFlorent Pruvost 80688c8ee7SFlorent Pruvost /* Use getRowIJ and the trick CSC/CSR instead of GetColumnIJ for performance */ 81688c8ee7SFlorent Pruvost PetscCall(MatGetRowIJ(A_loc, 0, PETSC_FALSE, PETSC_FALSE, &n, &row, &col, &flag)); 82688c8ee7SFlorent Pruvost PetscCheck(flag, PETSC_COMM_SELF, PETSC_ERR_SUP, "GetRowIJ failed"); 83688c8ee7SFlorent Pruvost PetscCall(MatSeqAIJGetArrayRead(A_loc, &val)); 84688c8ee7SFlorent Pruvost 85688c8ee7SFlorent Pruvost PetscCall(PetscMalloc1(1, &spm)); 86688c8ee7SFlorent Pruvost PetscStackCallExternalVoid("spmInitDist", spmInitDist(spm, pastix->comm)); 87688c8ee7SFlorent Pruvost 88688c8ee7SFlorent Pruvost spm->n = n; 89688c8ee7SFlorent Pruvost spm->nnz = row[n]; 90688c8ee7SFlorent Pruvost spm->fmttype = SpmCSR; 91688c8ee7SFlorent Pruvost spm->flttype = SPM_FLTTYPE; 92688c8ee7SFlorent Pruvost spm->replicated = !(A->rmap->n != A->rmap->N); 93688c8ee7SFlorent Pruvost 94688c8ee7SFlorent Pruvost PetscStackCallExternalVoid("spmUpdateComputedFields", spmUpdateComputedFields(spm)); 95688c8ee7SFlorent Pruvost PetscStackCallExternalVoid("spmAlloc", spmAlloc(spm)); 96688c8ee7SFlorent Pruvost 97688c8ee7SFlorent Pruvost /* Get data distribution */ 98688c8ee7SFlorent Pruvost if (!spm->replicated) { 99688c8ee7SFlorent Pruvost for (i = A->rmap->rstart; i < A->rmap->rend; i++) spm->loc2glob[i - A->rmap->rstart] = i; 100688c8ee7SFlorent Pruvost } 101688c8ee7SFlorent Pruvost 102688c8ee7SFlorent Pruvost /* Copy arrays */ 103688c8ee7SFlorent Pruvost PetscCall(PetscArraycpy(spm->colptr, col, spm->nnz)); 104688c8ee7SFlorent Pruvost PetscCall(PetscArraycpy(spm->rowptr, row, spm->n + 1)); 105688c8ee7SFlorent Pruvost PetscCall(PetscArraycpy((PetscScalar *)spm->values, val, spm->nnzexp)); 106688c8ee7SFlorent Pruvost PetscCall(MatRestoreRowIJ(A_loc, 0, PETSC_FALSE, PETSC_FALSE, &n, &row, &col, &flag)); 107688c8ee7SFlorent Pruvost PetscCheck(flag, PETSC_COMM_SELF, PETSC_ERR_SUP, "RestoreRowIJ failed"); 108688c8ee7SFlorent Pruvost PetscCall(MatSeqAIJRestoreArrayRead(A_loc, &val)); 109688c8ee7SFlorent Pruvost if (ismpiaij || ismpisbaij) PetscCall(MatDestroy(&A_loc)); 110688c8ee7SFlorent Pruvost 111688c8ee7SFlorent Pruvost /* 112688c8ee7SFlorent Pruvost PaStiX works only with CSC matrices for now, so let's lure him by telling him 113688c8ee7SFlorent Pruvost that the PETSc CSR matrix is a CSC matrix. 114688c8ee7SFlorent Pruvost Note that this is not available yet for Hermitian matrices and LL^h/LDL^h factorizations. 115745c78f7SBarry Smith */ 116688c8ee7SFlorent Pruvost { 117688c8ee7SFlorent Pruvost spm_int_t *tmp; 118688c8ee7SFlorent Pruvost spm->fmttype = SpmCSC; 119688c8ee7SFlorent Pruvost tmp = spm->colptr; 120688c8ee7SFlorent Pruvost spm->colptr = spm->rowptr; 121688c8ee7SFlorent Pruvost spm->rowptr = tmp; 122688c8ee7SFlorent Pruvost pastix->iparm[IPARM_TRANSPOSE_SOLVE] = PastixTrans; 123f31ce8a6SBarry Smith } 124745c78f7SBarry Smith 125688c8ee7SFlorent Pruvost /* Update matrix to be in PaStiX format */ 126688c8ee7SFlorent Pruvost PetscStackCallExternalVoid("spmCheckAndCorrect", spm_err = spmCheckAndCorrect(spm, &spm2)); 127688c8ee7SFlorent Pruvost if (spm_err != 0) { 128688c8ee7SFlorent Pruvost PetscStackCallExternalVoid("spmExit", spmExit(spm)); 129688c8ee7SFlorent Pruvost *spm = spm2; 1303bf14a46SMatthew Knepley } 1313bf14a46SMatthew Knepley 132688c8ee7SFlorent Pruvost if (isseqsbaij || ismpisbaij) PetscCall(MatDestroy(&A_aij)); 133688c8ee7SFlorent Pruvost 134688c8ee7SFlorent Pruvost pastix->spm = spm; 1353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1363bf14a46SMatthew Knepley } 1373bf14a46SMatthew Knepley 1383bf14a46SMatthew Knepley /* 139688c8ee7SFlorent Pruvost Call clean step of PaStiX if initialized 140688c8ee7SFlorent Pruvost Free the SpM matrix and the PaStiX instance. 1413bf14a46SMatthew Knepley */ 142688c8ee7SFlorent Pruvost static PetscErrorCode MatDestroy_PaStiX(Mat A) 143d71ae5a4SJacob Faibussowitsch { 144688c8ee7SFlorent Pruvost Mat_Pastix *pastix = (Mat_Pastix *)A->data; 145745c78f7SBarry Smith 1463bf14a46SMatthew Knepley PetscFunctionBegin; 147688c8ee7SFlorent Pruvost /* Finalize SPM (matrix handler of PaStiX) */ 148688c8ee7SFlorent Pruvost if (pastix->spm) { 149688c8ee7SFlorent Pruvost PetscStackCallExternalVoid("spmExit", spmExit(pastix->spm)); 150688c8ee7SFlorent Pruvost PetscCall(PetscFree(pastix->spm)); 1513bf14a46SMatthew Knepley } 152688c8ee7SFlorent Pruvost 153688c8ee7SFlorent Pruvost /* clear composed functions */ 1542e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorGetSolverType_C", NULL)); 155688c8ee7SFlorent Pruvost 156688c8ee7SFlorent Pruvost /* Finalize PaStiX */ 157688c8ee7SFlorent Pruvost if (pastix->pastix_data) pastixFinalize(&pastix->pastix_data); 158688c8ee7SFlorent Pruvost 159688c8ee7SFlorent Pruvost /* Deallocate PaStiX structure */ 1609566063dSJacob Faibussowitsch PetscCall(PetscFree(A->data)); 1613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1623bf14a46SMatthew Knepley } 1633bf14a46SMatthew Knepley 1643bf14a46SMatthew Knepley /* 165dd8e379bSPierre Jolivet Gather right-hand side. 1663bf14a46SMatthew Knepley Call for Solve step. 1673bf14a46SMatthew Knepley Scatter solution. 1683bf14a46SMatthew Knepley */ 16966976f2fSJacob Faibussowitsch static PetscErrorCode MatSolve_PaStiX(Mat A, Vec b, Vec x) 170d71ae5a4SJacob Faibussowitsch { 171688c8ee7SFlorent Pruvost Mat_Pastix *pastix = (Mat_Pastix *)A->data; 172688c8ee7SFlorent Pruvost const PetscScalar *bptr; 173688c8ee7SFlorent Pruvost PetscInt ldrhs; 1743bf14a46SMatthew Knepley 1753bf14a46SMatthew Knepley PetscFunctionBegin; 176688c8ee7SFlorent Pruvost pastix->rhsnbr = 1; 177688c8ee7SFlorent Pruvost ldrhs = pastix->spm->n; 178688c8ee7SFlorent Pruvost 1799566063dSJacob Faibussowitsch PetscCall(VecCopy(b, x)); 180688c8ee7SFlorent Pruvost PetscCall(VecGetArray(x, &pastix->rhs)); 181688c8ee7SFlorent Pruvost PetscCall(VecGetArrayRead(b, &bptr)); 1823bf14a46SMatthew Knepley 1833bf14a46SMatthew Knepley /* solve phase */ 184688c8ee7SFlorent Pruvost /*-------------*/ 185688c8ee7SFlorent Pruvost PetscCheck(pastix->pastix_data, PETSC_COMM_SELF, PETSC_ERR_SUP, "PaStiX hasn't been initialized"); 186688c8ee7SFlorent Pruvost PetscCallExternal(pastix_task_solve, pastix->pastix_data, ldrhs, pastix->rhsnbr, pastix->rhs, ldrhs); 187688c8ee7SFlorent Pruvost PetscCallExternal(pastix_task_refine, pastix->pastix_data, ldrhs, pastix->rhsnbr, (PetscScalar *)bptr, ldrhs, pastix->rhs, ldrhs); 1883bf14a46SMatthew Knepley 189688c8ee7SFlorent Pruvost PetscCall(VecRestoreArray(x, &pastix->rhs)); 190688c8ee7SFlorent Pruvost PetscCall(VecRestoreArrayRead(b, &bptr)); 1913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1923bf14a46SMatthew Knepley } 1933bf14a46SMatthew Knepley 1943bf14a46SMatthew Knepley /* 1953bf14a46SMatthew Knepley Numeric factorisation using PaStiX solver. 1963bf14a46SMatthew Knepley 197688c8ee7SFlorent Pruvost input: 198688c8ee7SFlorent Pruvost F - PETSc matrix that contains PaStiX interface. 199688c8ee7SFlorent Pruvost A - PETSC matrix in aij, bail or sbaij format 2003bf14a46SMatthew Knepley */ 20166976f2fSJacob Faibussowitsch static PetscErrorCode MatFactorNumeric_PaStiX(Mat F, Mat A, const MatFactorInfo *info) 202d71ae5a4SJacob Faibussowitsch { 203688c8ee7SFlorent Pruvost Mat_Pastix *pastix = (Mat_Pastix *)F->data; 2043bf14a46SMatthew Knepley 2053bf14a46SMatthew Knepley PetscFunctionBegin; 20657508eceSPierre Jolivet F->ops->solve = MatSolve_PaStiX; 2073bf14a46SMatthew Knepley 208688c8ee7SFlorent Pruvost /* Perform Numerical Factorization */ 209688c8ee7SFlorent Pruvost PetscCheck(pastix->pastix_data, PETSC_COMM_SELF, PETSC_ERR_SUP, "PaStiX hasn't been initialized"); 210688c8ee7SFlorent Pruvost PetscCallExternal(pastix_task_numfact, pastix->pastix_data, pastix->spm); 2113bf14a46SMatthew Knepley 21257508eceSPierre Jolivet F->assembled = PETSC_TRUE; 2133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2143bf14a46SMatthew Knepley } 2153bf14a46SMatthew Knepley 216688c8ee7SFlorent Pruvost static PetscErrorCode MatLUFactorNumeric_PaStiX(Mat F, Mat A, const MatFactorInfo *info) 217d71ae5a4SJacob Faibussowitsch { 218688c8ee7SFlorent Pruvost Mat_Pastix *pastix = (Mat_Pastix *)F->data; 2193bf14a46SMatthew Knepley 2203bf14a46SMatthew Knepley PetscFunctionBegin; 221688c8ee7SFlorent Pruvost PetscCheck(pastix->iparm[IPARM_FACTORIZATION] == PastixFactGETRF, PetscObjectComm((PetscObject)F), PETSC_ERR_SUP, "Incorrect factorization type for symbolic and numerical factorization by PaStiX"); 222688c8ee7SFlorent Pruvost pastix->iparm[IPARM_FACTORIZATION] = PastixFactGETRF; 223688c8ee7SFlorent Pruvost PetscCall(MatFactorNumeric_PaStiX(F, A, info)); 2243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2253bf14a46SMatthew Knepley } 2263bf14a46SMatthew Knepley 227688c8ee7SFlorent Pruvost static PetscErrorCode MatCholeskyFactorNumeric_PaStiX(Mat F, Mat A, const MatFactorInfo *info) 228d71ae5a4SJacob Faibussowitsch { 229688c8ee7SFlorent Pruvost Mat_Pastix *pastix = (Mat_Pastix *)F->data; 2303bf14a46SMatthew Knepley 2313bf14a46SMatthew Knepley PetscFunctionBegin; 232688c8ee7SFlorent Pruvost PetscCheck(pastix->iparm[IPARM_FACTORIZATION] == PastixFactSYTRF, PetscObjectComm((PetscObject)F), PETSC_ERR_SUP, "Incorrect factorization type for symbolic and numerical factorization by PaStiX"); 233688c8ee7SFlorent Pruvost pastix->iparm[IPARM_FACTORIZATION] = PastixFactSYTRF; 234688c8ee7SFlorent Pruvost PetscCall(MatFactorNumeric_PaStiX(F, A, info)); 235688c8ee7SFlorent Pruvost PetscFunctionReturn(PETSC_SUCCESS); 236688c8ee7SFlorent Pruvost } 237688c8ee7SFlorent Pruvost 238688c8ee7SFlorent Pruvost /* 239688c8ee7SFlorent Pruvost Perform Ordering step and Symbolic Factorization step 240688c8ee7SFlorent Pruvost 241688c8ee7SFlorent Pruvost Note the Petsc r and c permutations are ignored 242688c8ee7SFlorent Pruvost input: 243688c8ee7SFlorent Pruvost F - PETSc matrix that contains PaStiX interface. 244688c8ee7SFlorent Pruvost A - matrix in aij, bail or sbaij format 245688c8ee7SFlorent Pruvost r - permutation ? 246688c8ee7SFlorent Pruvost c - TODO 247*bfe80ac4SPierre Jolivet info - Information about the factorization to perform. 248688c8ee7SFlorent Pruvost output: 249688c8ee7SFlorent Pruvost pastix_data - This instance will be updated with the SolverMatrix allocated. 250688c8ee7SFlorent Pruvost */ 251688c8ee7SFlorent Pruvost static PetscErrorCode MatFactorSymbolic_PaStiX(Mat F, Mat A, IS r, IS c, const MatFactorInfo *info) 252688c8ee7SFlorent Pruvost { 253688c8ee7SFlorent Pruvost Mat_Pastix *pastix = (Mat_Pastix *)F->data; 254688c8ee7SFlorent Pruvost 255688c8ee7SFlorent Pruvost PetscFunctionBegin; 256688c8ee7SFlorent Pruvost pastix->matstruc = DIFFERENT_NONZERO_PATTERN; 257688c8ee7SFlorent Pruvost 258688c8ee7SFlorent Pruvost /* Initialise SPM structure */ 259688c8ee7SFlorent Pruvost PetscCall(MatConvertToSPM(A, MAT_INITIAL_MATRIX, pastix)); 260688c8ee7SFlorent Pruvost 261688c8ee7SFlorent Pruvost /* Ordering - Symbolic factorization - Build SolverMatrix */ 262688c8ee7SFlorent Pruvost PetscCheck(pastix->pastix_data, PETSC_COMM_SELF, PETSC_ERR_SUP, "PaStiX hasn't been initialized"); 263688c8ee7SFlorent Pruvost PetscCallExternal(pastix_task_analyze, pastix->pastix_data, pastix->spm); 264688c8ee7SFlorent Pruvost PetscFunctionReturn(PETSC_SUCCESS); 265688c8ee7SFlorent Pruvost } 266688c8ee7SFlorent Pruvost 267688c8ee7SFlorent Pruvost static PetscErrorCode MatLUFactorSymbolic_PaStiX(Mat F, Mat A, IS r, IS c, const MatFactorInfo *info) 268688c8ee7SFlorent Pruvost { 269688c8ee7SFlorent Pruvost Mat_Pastix *pastix = (Mat_Pastix *)F->data; 270688c8ee7SFlorent Pruvost 271688c8ee7SFlorent Pruvost PetscFunctionBegin; 272688c8ee7SFlorent Pruvost pastix->iparm[IPARM_FACTORIZATION] = PastixFactGETRF; 273688c8ee7SFlorent Pruvost PetscCall(MatFactorSymbolic_PaStiX(F, A, r, c, info)); 274688c8ee7SFlorent Pruvost PetscFunctionReturn(PETSC_SUCCESS); 275688c8ee7SFlorent Pruvost } 276688c8ee7SFlorent Pruvost 277688c8ee7SFlorent Pruvost /* Note the Petsc r permutation is ignored */ 278688c8ee7SFlorent Pruvost static PetscErrorCode MatCholeskyFactorSymbolic_PaStiX(Mat F, Mat A, IS r, const MatFactorInfo *info) 279688c8ee7SFlorent Pruvost { 280688c8ee7SFlorent Pruvost Mat_Pastix *pastix = (Mat_Pastix *)F->data; 281688c8ee7SFlorent Pruvost 282688c8ee7SFlorent Pruvost PetscFunctionBegin; 283688c8ee7SFlorent Pruvost /* Warning: Cholesky in Petsc wrapper does not handle (complex) Hermitian matrices. 284688c8ee7SFlorent Pruvost The factorization type can be forced using the parameter 285688c8ee7SFlorent Pruvost mat_pastix_factorization (see enum pastix_factotype_e in 286688c8ee7SFlorent Pruvost https://solverstack.gitlabpages.inria.fr/pastix/group__pastix__api.html). */ 287688c8ee7SFlorent Pruvost pastix->iparm[IPARM_FACTORIZATION] = PastixFactSYTRF; 288688c8ee7SFlorent Pruvost PetscCall(MatFactorSymbolic_PaStiX(F, A, r, NULL, info)); 2893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2903bf14a46SMatthew Knepley } 2913bf14a46SMatthew Knepley 29266976f2fSJacob Faibussowitsch static PetscErrorCode MatView_PaStiX(Mat A, PetscViewer viewer) 293d71ae5a4SJacob Faibussowitsch { 294ace3abfcSBarry Smith PetscBool iascii; 2953bf14a46SMatthew Knepley PetscViewerFormat format; 2963bf14a46SMatthew Knepley 2973bf14a46SMatthew Knepley PetscFunctionBegin; 2989566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 2993bf14a46SMatthew Knepley if (iascii) { 3009566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 3013bf14a46SMatthew Knepley if (format == PETSC_VIEWER_ASCII_INFO) { 302688c8ee7SFlorent Pruvost Mat_Pastix *pastix = (Mat_Pastix *)A->data; 303688c8ee7SFlorent Pruvost spmatrix_t *spm = pastix->spm; 304688c8ee7SFlorent Pruvost PetscCheck(!spm, PETSC_COMM_SELF, PETSC_ERR_SUP, "Sparse matrix isn't initialized"); 305b5e56a35SBarry Smith 3069566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "PaStiX run parameters:\n")); 307688c8ee7SFlorent Pruvost PetscCall(PetscViewerASCIIPrintf(viewer, " Matrix type : %s \n", ((spm->mtxtype == SpmSymmetric) ? "Symmetric" : "Unsymmetric"))); 308688c8ee7SFlorent Pruvost PetscCall(PetscViewerASCIIPrintf(viewer, " Level of printing (0,1,2): %ld \n", (long)pastix->iparm[IPARM_VERBOSE])); 309688c8ee7SFlorent Pruvost PetscCall(PetscViewerASCIIPrintf(viewer, " Number of refinements iterations : %ld \n", (long)pastix->iparm[IPARM_NBITER])); 310688c8ee7SFlorent Pruvost PetscCall(PetscPrintf(PETSC_COMM_SELF, " Error : %e \n", pastix->dparm[DPARM_RELATIVE_ERROR])); 311688c8ee7SFlorent Pruvost if (pastix->iparm[IPARM_VERBOSE] > 0) spmPrintInfo(spm, stdout); 3123bf14a46SMatthew Knepley } 3133bf14a46SMatthew Knepley } 3143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3153bf14a46SMatthew Knepley } 3163bf14a46SMatthew Knepley 3173bf14a46SMatthew Knepley /*MC 3182692d6eeSBarry Smith MATSOLVERPASTIX - A solver package providing direct solvers (LU) for distributed 3193bf14a46SMatthew Knepley and sequential matrices via the external package PaStiX. 3203bf14a46SMatthew Knepley 321688c8ee7SFlorent Pruvost Use `./configure --download-hwloc --download-metis --download-ptscotch --download-pastix --download-netlib-lapack [or MKL for ex. --with-blaslapack-dir=${MKLROOT}]` 322688c8ee7SFlorent Pruvost to have PETSc installed with PasTiX. 323c2b89b5dSBarry Smith 324688c8ee7SFlorent Pruvost Use `-pc_type lu` `-pc_factor_mat_solver_type pastix` to use this direct solver. 3253bf14a46SMatthew Knepley 3263bf14a46SMatthew Knepley Options Database Keys: 327688c8ee7SFlorent Pruvost -mat_pastix_verbose <0,1,2> - print level of information messages from PaStiX 328688c8ee7SFlorent Pruvost -mat_pastix_factorization <0,1,2,3,4> - Factorization algorithm (Cholesky, LDL^t, LU, LL^t, LDL^h) 329688c8ee7SFlorent Pruvost -mat_pastix_itermax <integer> - Maximum number of iterations during refinement 330688c8ee7SFlorent Pruvost -mat_pastix_epsilon_refinement <double> - Epsilon for refinement 331688c8ee7SFlorent Pruvost -mat_pastix_epsilon_magn_ctrl <double> - Epsilon for magnitude control 332688c8ee7SFlorent Pruvost -mat_pastix_ordering <0,1> - Ordering (Scotch or Metis) 333688c8ee7SFlorent Pruvost -mat_pastix_thread_nbr <integer> - Set the numbers of threads for each MPI process 334688c8ee7SFlorent Pruvost -mat_pastix_scheduler <0,1,2,3,4> - Scheduler (sequential, static, parsec, starpu, dynamic) 335688c8ee7SFlorent Pruvost -mat_pastix_compress_when <0,1,2,3> - When to compress (never, minimal-theory, just-in-time, supernodes) 336688c8ee7SFlorent Pruvost -mat_pastix_compress_tolerance <double> - Tolerance for low-rank kernels. 3373bf14a46SMatthew Knepley 33895452b02SPatrick Sanan Notes: 33995452b02SPatrick Sanan This only works for matrices with symmetric nonzero structure, if you pass it a matrix with 3402ef1f0ffSBarry Smith nonsymmetric structure PasTiX, and hence, PETSc return with an error. 341baed9decSBarry Smith 3423bf14a46SMatthew Knepley Level: beginner 3433bf14a46SMatthew Knepley 3441cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `PCFactorSetMatSolverType()`, `MatSolverType`, `MatGetFactor()` 3453bf14a46SMatthew Knepley M*/ 3463bf14a46SMatthew Knepley 34766976f2fSJacob Faibussowitsch static PetscErrorCode MatGetInfo_PaStiX(Mat A, MatInfoType flag, MatInfo *info) 348d71ae5a4SJacob Faibussowitsch { 349688c8ee7SFlorent Pruvost Mat_Pastix *pastix = (Mat_Pastix *)A->data; 3503bf14a46SMatthew Knepley 3513bf14a46SMatthew Knepley PetscFunctionBegin; 3523bf14a46SMatthew Knepley info->block_size = 1.0; 353688c8ee7SFlorent Pruvost info->nz_allocated = pastix->iparm[IPARM_ALLOCATED_TERMS]; 354688c8ee7SFlorent Pruvost info->nz_used = pastix->iparm[IPARM_NNZEROS]; 3553bf14a46SMatthew Knepley info->nz_unneeded = 0.0; 3563bf14a46SMatthew Knepley info->assemblies = 0.0; 3573bf14a46SMatthew Knepley info->mallocs = 0.0; 3583bf14a46SMatthew Knepley info->memory = 0.0; 3593bf14a46SMatthew Knepley info->fill_ratio_given = 0; 3603bf14a46SMatthew Knepley info->fill_ratio_needed = 0; 3613bf14a46SMatthew Knepley info->factor_mallocs = 0; 3623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3633bf14a46SMatthew Knepley } 3643bf14a46SMatthew Knepley 365688c8ee7SFlorent Pruvost static PetscErrorCode MatFactorGetSolverType_PaStiX(Mat A, MatSolverType *type) 366d71ae5a4SJacob Faibussowitsch { 3673bf14a46SMatthew Knepley PetscFunctionBegin; 3682692d6eeSBarry Smith *type = MATSOLVERPASTIX; 3693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3703bf14a46SMatthew Knepley } 3713bf14a46SMatthew Knepley 372688c8ee7SFlorent Pruvost /* Sets PaStiX options from the options database */ 373688c8ee7SFlorent Pruvost static PetscErrorCode MatSetFromOptions_PaStiX(Mat A) 374688c8ee7SFlorent Pruvost { 375688c8ee7SFlorent Pruvost Mat_Pastix *pastix = (Mat_Pastix *)A->data; 376688c8ee7SFlorent Pruvost pastix_int_t *iparm = pastix->iparm; 377688c8ee7SFlorent Pruvost double *dparm = pastix->dparm; 378688c8ee7SFlorent Pruvost PetscInt icntl; 379688c8ee7SFlorent Pruvost PetscReal dcntl; 380688c8ee7SFlorent Pruvost PetscBool set; 381688c8ee7SFlorent Pruvost 382688c8ee7SFlorent Pruvost PetscFunctionBegin; 383688c8ee7SFlorent Pruvost PetscOptionsBegin(PetscObjectComm((PetscObject)A), ((PetscObject)A)->prefix, "PaStiX Options", "Mat"); 384688c8ee7SFlorent Pruvost 385688c8ee7SFlorent Pruvost iparm[IPARM_VERBOSE] = 0; 386688c8ee7SFlorent Pruvost iparm[IPARM_ITERMAX] = 20; 387688c8ee7SFlorent Pruvost 388688c8ee7SFlorent Pruvost PetscCall(PetscOptionsRangeInt("-mat_pastix_verbose", "iparm[IPARM_VERBOSE] : level of printing (0 to 2)", "None", iparm[IPARM_VERBOSE], &icntl, &set, 0, 2)); 389688c8ee7SFlorent Pruvost if (set) iparm[IPARM_VERBOSE] = (pastix_int_t)icntl; 390688c8ee7SFlorent Pruvost 391688c8ee7SFlorent Pruvost PetscCall(PetscOptionsRangeInt("-mat_pastix_factorization", "iparm[IPARM_FACTORIZATION]: Factorization algorithm", "None", iparm[IPARM_FACTORIZATION], &icntl, &set, 0, 4)); 392688c8ee7SFlorent Pruvost if (set) iparm[IPARM_FACTORIZATION] = (pastix_int_t)icntl; 393688c8ee7SFlorent Pruvost 394688c8ee7SFlorent Pruvost PetscCall(PetscOptionsBoundedInt("-mat_pastix_itermax", "iparm[IPARM_ITERMAX]: Max iterations", "None", iparm[IPARM_ITERMAX], &icntl, &set, 1)); 395688c8ee7SFlorent Pruvost if (set) iparm[IPARM_ITERMAX] = (pastix_int_t)icntl; 396688c8ee7SFlorent Pruvost 397688c8ee7SFlorent Pruvost PetscCall(PetscOptionsBoundedReal("-mat_pastix_epsilon_refinement", "dparm[DPARM_EPSILON_REFINEMENT]: Epsilon refinement", "None", dparm[DPARM_EPSILON_REFINEMENT], &dcntl, &set, -1.)); 398688c8ee7SFlorent Pruvost if (set) dparm[DPARM_EPSILON_REFINEMENT] = (double)dcntl; 399688c8ee7SFlorent Pruvost 400688c8ee7SFlorent Pruvost PetscCall(PetscOptionsReal("-mat_pastix_epsilon_magn_ctrl", "dparm[DPARM_EPSILON_MAGN_CTRL]: Epsilon magnitude control", "None", dparm[DPARM_EPSILON_MAGN_CTRL], &dcntl, &set)); 401688c8ee7SFlorent Pruvost if (set) dparm[DPARM_EPSILON_MAGN_CTRL] = (double)dcntl; 402688c8ee7SFlorent Pruvost 403688c8ee7SFlorent Pruvost PetscCall(PetscOptionsRangeInt("-mat_pastix_ordering", "iparm[IPARM_ORDERING]: Ordering algorithm", "None", iparm[IPARM_ORDERING], &icntl, &set, 0, 2)); 404688c8ee7SFlorent Pruvost if (set) iparm[IPARM_ORDERING] = (pastix_int_t)icntl; 405688c8ee7SFlorent Pruvost 406688c8ee7SFlorent Pruvost PetscCall(PetscOptionsBoundedInt("-mat_pastix_thread_nbr", "iparm[IPARM_THREAD_NBR]: Number of thread by MPI node", "None", iparm[IPARM_THREAD_NBR], &icntl, &set, -1)); 407688c8ee7SFlorent Pruvost if (set) iparm[IPARM_THREAD_NBR] = (pastix_int_t)icntl; 408688c8ee7SFlorent Pruvost 409688c8ee7SFlorent Pruvost PetscCall(PetscOptionsRangeInt("-mat_pastix_scheduler", "iparm[IPARM_SCHEDULER]: Scheduler", "None", iparm[IPARM_SCHEDULER], &icntl, &set, 0, 4)); 410688c8ee7SFlorent Pruvost if (set) iparm[IPARM_SCHEDULER] = (pastix_int_t)icntl; 411688c8ee7SFlorent Pruvost 412688c8ee7SFlorent Pruvost PetscCall(PetscOptionsRangeInt("-mat_pastix_compress_when", "iparm[IPARM_COMPRESS_WHEN]: When to compress", "None", iparm[IPARM_COMPRESS_WHEN], &icntl, &set, 0, 3)); 413688c8ee7SFlorent Pruvost if (set) iparm[IPARM_COMPRESS_WHEN] = (pastix_int_t)icntl; 414688c8ee7SFlorent Pruvost 415688c8ee7SFlorent Pruvost PetscCall(PetscOptionsBoundedReal("-mat_pastix_compress_tolerance", "dparm[DPARM_COMPRESS_TOLERANCE]: Tolerance for low-rank kernels", "None", dparm[DPARM_COMPRESS_TOLERANCE], &dcntl, &set, 0.)); 416688c8ee7SFlorent Pruvost if (set) dparm[DPARM_COMPRESS_TOLERANCE] = (double)dcntl; 417688c8ee7SFlorent Pruvost 418688c8ee7SFlorent Pruvost PetscOptionsEnd(); 419688c8ee7SFlorent Pruvost PetscFunctionReturn(PETSC_SUCCESS); 420688c8ee7SFlorent Pruvost } 421688c8ee7SFlorent Pruvost 422688c8ee7SFlorent Pruvost static PetscErrorCode MatGetFactor_pastix(Mat A, MatFactorType ftype, Mat *F, const char *mattype) 423d71ae5a4SJacob Faibussowitsch { 4243bf14a46SMatthew Knepley Mat B; 4253bf14a46SMatthew Knepley Mat_Pastix *pastix; 4263bf14a46SMatthew Knepley 4273bf14a46SMatthew Knepley PetscFunctionBegin; 4283bf14a46SMatthew Knepley /* Create the factorization matrix */ 4299566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B)); 4309566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N)); 431688c8ee7SFlorent Pruvost PetscCall(PetscStrallocpy(MATSOLVERPASTIX, &((PetscObject)B)->type_name)); 4329566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 4333bf14a46SMatthew Knepley 434688c8ee7SFlorent Pruvost PetscCheck(ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_CHOLESKY, PETSC_COMM_SELF, PETSC_ERR_SUP, "Factor type not supported by PaStiX"); 4353bf14a46SMatthew Knepley 43600c67f3bSHong Zhang /* set solvertype */ 4379566063dSJacob Faibussowitsch PetscCall(PetscFree(B->solvertype)); 4389566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSOLVERPASTIX, &B->solvertype)); 43900c67f3bSHong Zhang 440688c8ee7SFlorent Pruvost B->ops->lufactorsymbolic = MatLUFactorSymbolic_PaStiX; 441688c8ee7SFlorent Pruvost B->ops->lufactornumeric = MatLUFactorNumeric_PaStiX; 442688c8ee7SFlorent Pruvost B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_PaStiX; 443688c8ee7SFlorent Pruvost B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_PaStiX; 444688c8ee7SFlorent Pruvost B->ops->view = MatView_PaStiX; 445688c8ee7SFlorent Pruvost B->ops->getinfo = MatGetInfo_PaStiX; 446688c8ee7SFlorent Pruvost B->ops->destroy = MatDestroy_PaStiX; 4472205254eSKarl Rupp 448688c8ee7SFlorent Pruvost PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_PaStiX)); 449688c8ee7SFlorent Pruvost 450688c8ee7SFlorent Pruvost B->factortype = ftype; 451688c8ee7SFlorent Pruvost 452688c8ee7SFlorent Pruvost /* Create the pastix structure */ 453688c8ee7SFlorent Pruvost PetscCall(PetscNew(&pastix)); 45458c95f1bSBarry Smith B->data = (void *)pastix; 4553bf14a46SMatthew Knepley 456688c8ee7SFlorent Pruvost /* Call to set default pastix options */ 457688c8ee7SFlorent Pruvost PetscStackCallExternalVoid("pastixInitParam", pastixInitParam(pastix->iparm, pastix->dparm)); 458688c8ee7SFlorent Pruvost PetscCall(MatSetFromOptions_PaStiX(B)); 459688c8ee7SFlorent Pruvost 460688c8ee7SFlorent Pruvost /* Get PETSc communicator */ 461688c8ee7SFlorent Pruvost PetscCall(PetscObjectGetComm((PetscObject)A, &pastix->comm)); 462688c8ee7SFlorent Pruvost 463688c8ee7SFlorent Pruvost /* Initialise PaStiX structure */ 464688c8ee7SFlorent Pruvost pastix->iparm[IPARM_SCOTCH_MT] = 0; 465688c8ee7SFlorent Pruvost PetscStackCallExternalVoid("pastixInit", pastixInit(&pastix->pastix_data, pastix->comm, pastix->iparm, pastix->dparm)); 466688c8ee7SFlorent Pruvost 467688c8ee7SFlorent Pruvost /* Warning: Cholesky in Petsc wrapper does not handle (complex) Hermitian matrices. 468688c8ee7SFlorent Pruvost The factorization type can be forced using the parameter 469688c8ee7SFlorent Pruvost mat_pastix_factorization (see enum pastix_factotype_e in 470688c8ee7SFlorent Pruvost https://solverstack.gitlabpages.inria.fr/pastix/group__pastix__api.html). */ 471688c8ee7SFlorent Pruvost pastix->iparm[IPARM_FACTORIZATION] = ftype == MAT_FACTOR_CHOLESKY ? PastixFactSYTRF : PastixFactGETRF; 472688c8ee7SFlorent Pruvost 4733bf14a46SMatthew Knepley *F = B; 4743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4753bf14a46SMatthew Knepley } 4763bf14a46SMatthew Knepley 477d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetFactor_mpiaij_pastix(Mat A, MatFactorType ftype, Mat *F) 478d71ae5a4SJacob Faibussowitsch { 4793bf14a46SMatthew Knepley PetscFunctionBegin; 4805f80ce2aSJacob Faibussowitsch PetscCheck(ftype == MAT_FACTOR_LU, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot use PETSc AIJ matrices with PaStiX Cholesky, use SBAIJ matrix"); 481688c8ee7SFlorent Pruvost PetscCall(MatGetFactor_pastix(A, ftype, F, MATMPIAIJ)); 4823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4833bf14a46SMatthew Knepley } 4843bf14a46SMatthew Knepley 485688c8ee7SFlorent Pruvost static PetscErrorCode MatGetFactor_seqaij_pastix(Mat A, MatFactorType ftype, Mat *F) 486d71ae5a4SJacob Faibussowitsch { 4873bf14a46SMatthew Knepley PetscFunctionBegin; 488688c8ee7SFlorent Pruvost PetscCheck(ftype == MAT_FACTOR_LU, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot use PETSc AIJ matrices with PaStiX Cholesky, use SBAIJ matrix"); 489688c8ee7SFlorent Pruvost PetscCall(MatGetFactor_pastix(A, ftype, F, MATSEQAIJ)); 4903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4913bf14a46SMatthew Knepley } 4923bf14a46SMatthew Knepley 493d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetFactor_mpisbaij_pastix(Mat A, MatFactorType ftype, Mat *F) 494d71ae5a4SJacob Faibussowitsch { 4953bf14a46SMatthew Knepley PetscFunctionBegin; 4965f80ce2aSJacob Faibussowitsch PetscCheck(ftype == MAT_FACTOR_CHOLESKY, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot use PETSc SBAIJ matrices with PaStiX LU, use AIJ matrix"); 497688c8ee7SFlorent Pruvost PetscCall(MatGetFactor_pastix(A, ftype, F, MATMPISBAIJ)); 498688c8ee7SFlorent Pruvost PetscFunctionReturn(PETSC_SUCCESS); 499688c8ee7SFlorent Pruvost } 50041c8de11SBarry Smith 501688c8ee7SFlorent Pruvost static PetscErrorCode MatGetFactor_seqsbaij_pastix(Mat A, MatFactorType ftype, Mat *F) 502688c8ee7SFlorent Pruvost { 503688c8ee7SFlorent Pruvost PetscFunctionBegin; 504688c8ee7SFlorent Pruvost PetscCheck(ftype == MAT_FACTOR_CHOLESKY, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot use PETSc SBAIJ matrices with PaStiX LU, use AIJ matrix"); 505688c8ee7SFlorent Pruvost PetscCall(MatGetFactor_pastix(A, ftype, F, MATSEQSBAIJ)); 5063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5073bf14a46SMatthew Knepley } 508f7a08781SBarry Smith 509d1f0640dSPierre Jolivet PETSC_INTERN PetscErrorCode MatSolverTypeRegister_Pastix(void) 510d71ae5a4SJacob Faibussowitsch { 51142c9c57cSBarry Smith PetscFunctionBegin; 5129566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERPASTIX, MATMPIAIJ, MAT_FACTOR_LU, MatGetFactor_mpiaij_pastix)); 5139566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERPASTIX, MATSEQAIJ, MAT_FACTOR_LU, MatGetFactor_seqaij_pastix)); 5149566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERPASTIX, MATMPISBAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_mpisbaij_pastix)); 5159566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERPASTIX, MATSEQSBAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_seqsbaij_pastix)); 5163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 51742c9c57cSBarry Smith } 518