148f88336SPieter Ghysels #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 27d6ea485SPieter Ghysels #include <../src/mat/impls/aij/mpi/mpiaij.h> 37d6ea485SPieter Ghysels #include <StrumpackSparseSolver.h> 47d6ea485SPieter Ghysels 5d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetDiagonal_STRUMPACK(Mat A, Vec v) 6d71ae5a4SJacob Faibussowitsch { 77d6ea485SPieter Ghysels PetscFunctionBegin; 87d6ea485SPieter Ghysels SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Mat type: STRUMPACK factor"); 93ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 107d6ea485SPieter Ghysels } 117d6ea485SPieter Ghysels 12d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_STRUMPACK(Mat A) 13d71ae5a4SJacob Faibussowitsch { 1435e8bcc9SJunchao Zhang STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)A->data; 157d6ea485SPieter Ghysels 167d6ea485SPieter Ghysels PetscFunctionBegin; 177d6ea485SPieter Ghysels /* Deallocate STRUMPACK storage */ 18e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_destroy", STRUMPACK_destroy(S)); 1935e8bcc9SJunchao Zhang PetscCall(PetscFree(A->data)); 20575704cbSPieter Ghysels 21575704cbSPieter Ghysels /* clear composed functions */ 229566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorGetSolverType_C", NULL)); 239566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetReordering_C", NULL)); 2429e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKGetReordering_C", NULL)); 259566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetColPerm_C", NULL)); 2629e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKGetColPerm_C", NULL)); 2729e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetGeometricNxyz_C", NULL)); 2829e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetGeometricComponents_C", NULL)); 2929e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetGeometricWidth_C", NULL)); 3029e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetGPU_C", NULL)); 3129e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKGetGPU_C", NULL)); 3229e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetCompression_C", NULL)); 3329e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKGetCompression_C", NULL)); 3429e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetCompRelTol_C", NULL)); 3529e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKGetCompRelTol_C", NULL)); 3629e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetCompAbsTol_C", NULL)); 3729e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKGetCompAbsTol_C", NULL)); 3829e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetCompMaxRank_C", NULL)); 3929e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKGetCompMaxRank_C", NULL)); 4029e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetCompLeafSize_C", NULL)); 4129e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKGetCompLeafSize_C", NULL)); 4229e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetCompMinSepSize_C", NULL)); 4329e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKGetCompMinSepSize_C", NULL)); 4429e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetCompLossyPrecision_C", NULL)); 4529e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKGetCompLossyPrecision_C", NULL)); 4629e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetCompButterflyLevels_C", NULL)); 4729e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKGetCompButterflyLevels_C", NULL)); 483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 49575704cbSPieter Ghysels } 50575704cbSPieter Ghysels 51d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSTRUMPACKSetReordering_STRUMPACK(Mat F, MatSTRUMPACKReordering reordering) 52d71ae5a4SJacob Faibussowitsch { 5335e8bcc9SJunchao Zhang STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 54575704cbSPieter Ghysels 55575704cbSPieter Ghysels PetscFunctionBegin; 5629e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_reordering_method", STRUMPACK_set_reordering_method(*S, (STRUMPACK_REORDERING_STRATEGY)reordering)); 5729e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 5829e0a805SPieter Ghysels } 5929e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetReordering_STRUMPACK(Mat F, MatSTRUMPACKReordering *reordering) 6029e0a805SPieter Ghysels { 6129e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 6229e0a805SPieter Ghysels 6329e0a805SPieter Ghysels PetscFunctionBegin; 6429e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_reordering_method", *reordering = (MatSTRUMPACKReordering)STRUMPACK_reordering_method(*S)); 653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6634f43fa5SPieter Ghysels } 67e5a36eccSStefano Zampini 68542aee0fSPieter Ghysels /*@ 691d27aa22SBarry Smith MatSTRUMPACKSetReordering - Set STRUMPACK <https://portal.nersc.gov/project/sparse/strumpack/master> fill-reducing reordering 7034f43fa5SPieter Ghysels 7129e0a805SPieter Ghysels Logically Collective 7229e0a805SPieter Ghysels 7334f43fa5SPieter Ghysels Input Parameters: 7429e0a805SPieter Ghysels + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 7529e0a805SPieter Ghysels - reordering - the code to be used to find the fill-reducing reordering 7634f43fa5SPieter Ghysels 7711a5261eSBarry Smith Options Database Key: 782ef1f0ffSBarry Smith . -mat_strumpack_reordering <METIS> - Sparsity reducing matrix reordering, see `MatSTRUMPACKReordering` 7934f43fa5SPieter Ghysels 8029e0a805SPieter Ghysels Level: intermediate 8134f43fa5SPieter Ghysels 821d27aa22SBarry Smith .seealso: `MATSOLVERSTRUMPACK`, [](ch_matrices), `Mat`, `MatSTRUMPACKReordering`, `MatGetFactor()`, `MatSTRUMPACKSetColPerm()`, `MatSTRUMPACKGetReordering()` 83542aee0fSPieter Ghysels @*/ 84d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSTRUMPACKSetReordering(Mat F, MatSTRUMPACKReordering reordering) 85d71ae5a4SJacob Faibussowitsch { 8634f43fa5SPieter Ghysels PetscFunctionBegin; 8734f43fa5SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 8834f43fa5SPieter Ghysels PetscValidLogicalCollectiveEnum(F, reordering, 2); 89cac4c232SBarry Smith PetscTryMethod(F, "MatSTRUMPACKSetReordering_C", (Mat, MatSTRUMPACKReordering), (F, reordering)); 903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 91575704cbSPieter Ghysels } 9229e0a805SPieter Ghysels /*@ 931d27aa22SBarry Smith MatSTRUMPACKGetReordering - Get STRUMPACK <https://portal.nersc.gov/project/sparse/strumpack/master> fill-reducing reordering 9429e0a805SPieter Ghysels 9529e0a805SPieter Ghysels Logically Collective 9629e0a805SPieter Ghysels 9729e0a805SPieter Ghysels Input Parameters: 9829e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 9929e0a805SPieter Ghysels 10029e0a805SPieter Ghysels Output Parameter: 10129e0a805SPieter Ghysels . reordering - the code to be used to find the fill-reducing reordering 10229e0a805SPieter Ghysels 10329e0a805SPieter Ghysels Level: intermediate 10429e0a805SPieter Ghysels 1051d27aa22SBarry Smith .seealso: `MATSOLVERSTRUMPACK`, [](ch_matrices), `Mat`, `MatSTRUMPACKReordering`, `MatGetFactor()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()` 10629e0a805SPieter Ghysels @*/ 10729e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetReordering(Mat F, MatSTRUMPACKReordering *reordering) 10829e0a805SPieter Ghysels { 10929e0a805SPieter Ghysels PetscFunctionBegin; 11029e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 11129e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKGetReordering_C", (Mat, MatSTRUMPACKReordering *), (F, reordering)); 11229e0a805SPieter Ghysels PetscValidLogicalCollectiveEnum(F, *reordering, 2); 11329e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 11429e0a805SPieter Ghysels } 115575704cbSPieter Ghysels 116d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSTRUMPACKSetColPerm_STRUMPACK(Mat F, PetscBool cperm) 117d71ae5a4SJacob Faibussowitsch { 11835e8bcc9SJunchao Zhang STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 11934f43fa5SPieter Ghysels 12034f43fa5SPieter Ghysels PetscFunctionBegin; 12129e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_matching", STRUMPACK_set_matching(*S, cperm ? STRUMPACK_MATCHING_MAX_DIAGONAL_PRODUCT_SCALING : STRUMPACK_MATCHING_NONE)); 12229e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 12329e0a805SPieter Ghysels } 12429e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetColPerm_STRUMPACK(Mat F, PetscBool *cperm) 12529e0a805SPieter Ghysels { 12629e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 12729e0a805SPieter Ghysels 12829e0a805SPieter Ghysels PetscFunctionBegin; 12929e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_matching", *cperm = (PetscBool)(STRUMPACK_matching(*S) != STRUMPACK_MATCHING_NONE)); 1303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13134f43fa5SPieter Ghysels } 132e5a36eccSStefano Zampini 133575704cbSPieter Ghysels /*@ 1341d27aa22SBarry Smith MatSTRUMPACKSetColPerm - Set whether STRUMPACK <https://portal.nersc.gov/project/sparse/strumpack/master> 1351d27aa22SBarry Smith should try to permute the columns of the matrix in order to get a nonzero diagonal 136147403d9SBarry Smith 137c3339decSBarry Smith Logically Collective 138575704cbSPieter Ghysels 139575704cbSPieter Ghysels Input Parameters: 1402ef1f0ffSBarry Smith + F - the factored matrix obtained by calling `MatGetFactor()` 14111a5261eSBarry Smith - cperm - `PETSC_TRUE` to permute (internally) the columns of the matrix 142575704cbSPieter Ghysels 14311a5261eSBarry Smith Options Database Key: 144147403d9SBarry Smith . -mat_strumpack_colperm <cperm> - true to use the permutation 145575704cbSPieter Ghysels 14629e0a805SPieter Ghysels Level: intermediate 147575704cbSPieter Ghysels 1481d27aa22SBarry Smith .seealso: `MATSOLVERSTRUMPACK`, [](ch_matrices), `MatSTRUMPACKSetReordering()`, `Mat`, `MatGetFactor()`, `MatSTRUMPACKGetColPerm()` 149575704cbSPieter Ghysels @*/ 150d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSTRUMPACKSetColPerm(Mat F, PetscBool cperm) 151d71ae5a4SJacob Faibussowitsch { 152575704cbSPieter Ghysels PetscFunctionBegin; 153575704cbSPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 15434f43fa5SPieter Ghysels PetscValidLogicalCollectiveBool(F, cperm, 2); 155cac4c232SBarry Smith PetscTryMethod(F, "MatSTRUMPACKSetColPerm_C", (Mat, PetscBool), (F, cperm)); 1563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 157575704cbSPieter Ghysels } 15829e0a805SPieter Ghysels /*@ 1591d27aa22SBarry Smith MatSTRUMPACKGetColPerm - Get whether STRUMPACK <https://portal.nersc.gov/project/sparse/strumpack/master> 1601d27aa22SBarry Smith will try to permute the columns of the matrix in order to get a nonzero diagonal 161575704cbSPieter Ghysels 16229e0a805SPieter Ghysels Logically Collective 16329e0a805SPieter Ghysels 16429e0a805SPieter Ghysels Input Parameters: 16529e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` 16629e0a805SPieter Ghysels 16729e0a805SPieter Ghysels Output Parameter: 16829e0a805SPieter Ghysels . cperm - Indicates whether STRUMPACK will permute columns 16929e0a805SPieter Ghysels 17029e0a805SPieter Ghysels Level: intermediate 17129e0a805SPieter Ghysels 1721d27aa22SBarry Smith .seealso: `MATSOLVERSTRUMPACK`, [](ch_matrices), `MatSTRUMPACKSetReordering()`, `Mat`, `MatGetFactor()`, `MatSTRUMPACKSetColPerm()` 17329e0a805SPieter Ghysels @*/ 17429e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetColPerm(Mat F, PetscBool *cperm) 17529e0a805SPieter Ghysels { 17629e0a805SPieter Ghysels PetscFunctionBegin; 17729e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 17829e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKGetColPerm_C", (Mat, PetscBool *), (F, cperm)); 17929e0a805SPieter Ghysels PetscValidLogicalCollectiveBool(F, *cperm, 2); 18029e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 18129e0a805SPieter Ghysels } 18229e0a805SPieter Ghysels 18329e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetGPU_STRUMPACK(Mat F, PetscBool gpu) 184d71ae5a4SJacob Faibussowitsch { 18535e8bcc9SJunchao Zhang STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 18634f43fa5SPieter Ghysels 18734f43fa5SPieter Ghysels PetscFunctionBegin; 18829e0a805SPieter Ghysels if (gpu) { 18929e0a805SPieter Ghysels #if !(defined(STRUMPACK_USE_CUDA) || defined(STRUMPACK_USE_HIP) || defined(STRUMPACK_USE_SYCL)) 19029e0a805SPieter Ghysels PetscCall(PetscPrintf(PETSC_COMM_SELF, " Warning: strumpack was not configured with GPU support\n")); 19129e0a805SPieter Ghysels #endif 19229e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_enable_gpu", STRUMPACK_enable_gpu(*S)); 19329e0a805SPieter Ghysels } else PetscStackCallExternalVoid("STRUMPACK_disable_gpu", STRUMPACK_disable_gpu(*S)); 19429e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 19529e0a805SPieter Ghysels } 19629e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetGPU_STRUMPACK(Mat F, PetscBool *gpu) 19729e0a805SPieter Ghysels { 19829e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 19929e0a805SPieter Ghysels 20029e0a805SPieter Ghysels PetscFunctionBegin; 20129e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_use_gpu", *gpu = (PetscBool)STRUMPACK_use_gpu(*S)); 2023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 20334f43fa5SPieter Ghysels } 204e5a36eccSStefano Zampini 20534f43fa5SPieter Ghysels /*@ 2061d27aa22SBarry Smith MatSTRUMPACKSetGPU - Set whether STRUMPACK <https://portal.nersc.gov/project/sparse/strumpack/master> 2071d27aa22SBarry Smith should enable GPU acceleration (not supported for all compression types) 20829e0a805SPieter Ghysels 20929e0a805SPieter Ghysels Logically Collective 21029e0a805SPieter Ghysels 21129e0a805SPieter Ghysels Input Parameters: 21229e0a805SPieter Ghysels + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 21329e0a805SPieter Ghysels - gpu - whether or not to use GPU acceleration 21429e0a805SPieter Ghysels 21529e0a805SPieter Ghysels Options Database Key: 21629e0a805SPieter Ghysels . -mat_strumpack_gpu <gpu> - true to use gpu offload 21729e0a805SPieter Ghysels 21829e0a805SPieter Ghysels Level: intermediate 21929e0a805SPieter Ghysels 2201d27aa22SBarry Smith .seealso: `MATSOLVERSTRUMPACK`, `MatGetFactor()`, `MatSTRUMPACKGetGPU()` 22129e0a805SPieter Ghysels @*/ 22229e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetGPU(Mat F, PetscBool gpu) 22329e0a805SPieter Ghysels { 22429e0a805SPieter Ghysels PetscFunctionBegin; 22529e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 22629e0a805SPieter Ghysels PetscValidLogicalCollectiveBool(F, gpu, 2); 22729e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKSetGPU_C", (Mat, PetscBool), (F, gpu)); 22829e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 22929e0a805SPieter Ghysels } 23029e0a805SPieter Ghysels /*@ 2311d27aa22SBarry Smith MatSTRUMPACKGetGPU - Get whether STRUMPACK <https://portal.nersc.gov/project/sparse/strumpack/master> 2321d27aa22SBarry Smith will try to use GPU acceleration (not supported for all compression types) 23329e0a805SPieter Ghysels 23429e0a805SPieter Ghysels Logically Collective 23529e0a805SPieter Ghysels 23629e0a805SPieter Ghysels Input Parameters: 23729e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 23829e0a805SPieter Ghysels 23929e0a805SPieter Ghysels Output Parameter: 24029e0a805SPieter Ghysels . gpu - whether or not STRUMPACK will try to use GPU acceleration 24129e0a805SPieter Ghysels 24229e0a805SPieter Ghysels Level: intermediate 24329e0a805SPieter Ghysels 2441d27aa22SBarry Smith .seealso: `MATSOLVERSTRUMPACK`, `MatGetFactor()`, `MatSTRUMPACKSetGPU()` 24529e0a805SPieter Ghysels @*/ 24629e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetGPU(Mat F, PetscBool *gpu) 24729e0a805SPieter Ghysels { 24829e0a805SPieter Ghysels PetscFunctionBegin; 24929e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 25029e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKGetGPU_C", (Mat, PetscBool *), (F, gpu)); 25129e0a805SPieter Ghysels PetscValidLogicalCollectiveBool(F, *gpu, 2); 25229e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 25329e0a805SPieter Ghysels } 25429e0a805SPieter Ghysels 25529e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetCompression_STRUMPACK(Mat F, MatSTRUMPACKCompressionType comp) 25629e0a805SPieter Ghysels { 25729e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 25829e0a805SPieter Ghysels 25929e0a805SPieter Ghysels PetscFunctionBegin; 260c2935fbeSBarry Smith #if !defined(STRUMPACK_USE_BPACK) 26129e0a805SPieter Ghysels PetscCheck(comp != MAT_STRUMPACK_COMPRESSION_TYPE_HODLR && comp != MAT_STRUMPACK_COMPRESSION_TYPE_BLR_HODLR && comp != MAT_STRUMPACK_COMPRESSION_TYPE_ZFP_BLR_HODLR, PetscObjectComm((PetscObject)F), PETSC_ERR_SUP, "Compression scheme requires ButterflyPACK, please reconfigure with --download-butterflypack"); 26229e0a805SPieter Ghysels #endif 263c2935fbeSBarry Smith #if !defined(STRUMPACK_USE_ZFP) 26429e0a805SPieter Ghysels PetscCheck(comp != MAT_STRUMPACK_COMPRESSION_TYPE_ZFP_BLR_HODLR && comp != MAT_STRUMPACK_COMPRESSION_TYPE_LOSSLESS && comp != MAT_STRUMPACK_COMPRESSION_TYPE_LOSSY, PetscObjectComm((PetscObject)F), PETSC_ERR_SUP, "Compression scheme requires ZFP, please reconfigure with --download-zfp"); 26529e0a805SPieter Ghysels #endif 26629e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_compression", STRUMPACK_set_compression(*S, (STRUMPACK_COMPRESSION_TYPE)comp)); 26729e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 26829e0a805SPieter Ghysels } 26929e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetCompression_STRUMPACK(Mat F, MatSTRUMPACKCompressionType *comp) 27029e0a805SPieter Ghysels { 27129e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 27229e0a805SPieter Ghysels 27329e0a805SPieter Ghysels PetscFunctionBegin; 27429e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression", *comp = (MatSTRUMPACKCompressionType)STRUMPACK_compression(*S)); 27529e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 27629e0a805SPieter Ghysels } 27729e0a805SPieter Ghysels 27829e0a805SPieter Ghysels /*@ 2791d27aa22SBarry Smith MatSTRUMPACKSetCompression - Set STRUMPACK <https://portal.nersc.gov/project/sparse/strumpack/master> compression type 28029e0a805SPieter Ghysels 28129e0a805SPieter Ghysels Input Parameters: 28229e0a805SPieter Ghysels + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 28329e0a805SPieter Ghysels - comp - Type of compression to be used in the approximate sparse factorization 28429e0a805SPieter Ghysels 28529e0a805SPieter Ghysels Options Database Key: 28629e0a805SPieter Ghysels . -mat_strumpack_compression <NONE> - Type of rank-structured compression in sparse LU factors (choose one of) NONE HSS BLR HODLR BLR_HODLR ZFP_BLR_HODLR LOSSLESS LOSSY 28729e0a805SPieter Ghysels 28829e0a805SPieter Ghysels Level: intermediate 28929e0a805SPieter Ghysels 2901d27aa22SBarry Smith Note: 29129e0a805SPieter Ghysels Default for `comp` is `MAT_STRUMPACK_COMPRESSION_TYPE_NONE` for `-pc_type lu` and `MAT_STRUMPACK_COMPRESSION_TYPE_BLR` 29229e0a805SPieter Ghysels for `-pc_type ilu` 29329e0a805SPieter Ghysels 2941d27aa22SBarry Smith .seealso: `MATSOLVERSTRUMPACK`, `MatGetFactor()`, `MatSTRUMPACKCompressionType`, `MatSTRUMPACKGetCompression()` 29529e0a805SPieter Ghysels @*/ 29629e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetCompression(Mat F, MatSTRUMPACKCompressionType comp) 29729e0a805SPieter Ghysels { 29829e0a805SPieter Ghysels PetscFunctionBegin; 29929e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 30029e0a805SPieter Ghysels PetscValidLogicalCollectiveEnum(F, comp, 2); 30129e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKSetCompression_C", (Mat, MatSTRUMPACKCompressionType), (F, comp)); 30229e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 30329e0a805SPieter Ghysels } 30429e0a805SPieter Ghysels /*@ 3051d27aa22SBarry Smith MatSTRUMPACKGetCompression - Get STRUMPACK <https://portal.nersc.gov/project/sparse/strumpack/master> compression type 30629e0a805SPieter Ghysels 30729e0a805SPieter Ghysels Input Parameters: 30829e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 30929e0a805SPieter Ghysels 31029e0a805SPieter Ghysels Output Parameter: 31129e0a805SPieter Ghysels . comp - Type of compression to be used in the approximate sparse factorization 31229e0a805SPieter Ghysels 31329e0a805SPieter Ghysels Level: intermediate 31429e0a805SPieter Ghysels 3151d27aa22SBarry Smith Note: 31629e0a805SPieter Ghysels Default is `MAT_STRUMPACK_COMPRESSION_TYPE_NONE` for `-pc_type lu` and `MAT_STRUMPACK_COMPRESSION_TYPE_BLR` for `-pc_type ilu` 31729e0a805SPieter Ghysels 3181d27aa22SBarry Smith .seealso: `MATSOLVERSTRUMPACK`, `MatGetFactor()`, `MatSTRUMPACKCompressionType`, `MatSTRUMPACKSetCompression()` 31929e0a805SPieter Ghysels @*/ 32029e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetCompression(Mat F, MatSTRUMPACKCompressionType *comp) 32129e0a805SPieter Ghysels { 32229e0a805SPieter Ghysels PetscFunctionBegin; 32329e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 32429e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKGetCompression_C", (Mat, MatSTRUMPACKCompressionType *), (F, comp)); 32529e0a805SPieter Ghysels PetscValidLogicalCollectiveEnum(F, *comp, 2); 32629e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 32729e0a805SPieter Ghysels } 32829e0a805SPieter Ghysels 32929e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetCompRelTol_STRUMPACK(Mat F, PetscReal rtol) 33029e0a805SPieter Ghysels { 33129e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 33229e0a805SPieter Ghysels 33329e0a805SPieter Ghysels PetscFunctionBegin; 33429e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_compression_rel_tol", STRUMPACK_set_compression_rel_tol(*S, rtol)); 33529e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 33629e0a805SPieter Ghysels } 33729e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetCompRelTol_STRUMPACK(Mat F, PetscReal *rtol) 33829e0a805SPieter Ghysels { 33929e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 34029e0a805SPieter Ghysels 34129e0a805SPieter Ghysels PetscFunctionBegin; 34229e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression_rel_tol", *rtol = (PetscReal)STRUMPACK_compression_rel_tol(*S)); 34329e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 34429e0a805SPieter Ghysels } 34529e0a805SPieter Ghysels 34629e0a805SPieter Ghysels /*@ 3471d27aa22SBarry Smith MatSTRUMPACKSetCompRelTol - Set STRUMPACK <https://portal.nersc.gov/project/sparse/strumpack/master> relative tolerance for compression 348147403d9SBarry Smith 349c3339decSBarry Smith Logically Collective 35034f43fa5SPieter Ghysels 35134f43fa5SPieter Ghysels Input Parameters: 3522ef1f0ffSBarry Smith + F - the factored matrix obtained by calling `MatGetFactor()` 35334f43fa5SPieter Ghysels - rtol - relative compression tolerance 35434f43fa5SPieter Ghysels 35511a5261eSBarry Smith Options Database Key: 35629e0a805SPieter Ghysels . -mat_strumpack_compression_rel_tol <1e-4> - Relative compression tolerance, when using `-pctype ilu` 35734f43fa5SPieter Ghysels 35829e0a805SPieter Ghysels Level: intermediate 35934f43fa5SPieter Ghysels 3601d27aa22SBarry Smith .seealso: `MATSOLVERSTRUMPACK`, [](ch_matrices), `Mat`, `MatGetFactor()`, `MatSTRUMPACKGetCompRelTol()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()` 36134f43fa5SPieter Ghysels @*/ 36229e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetCompRelTol(Mat F, PetscReal rtol) 363d71ae5a4SJacob Faibussowitsch { 36434f43fa5SPieter Ghysels PetscFunctionBegin; 36534f43fa5SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 36634f43fa5SPieter Ghysels PetscValidLogicalCollectiveReal(F, rtol, 2); 36729e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKSetCompRelTol_C", (Mat, PetscReal), (F, rtol)); 36829e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 36929e0a805SPieter Ghysels } 37029e0a805SPieter Ghysels /*@ 3711d27aa22SBarry Smith MatSTRUMPACKGetCompRelTol - Get STRUMPACK <https://portal.nersc.gov/project/sparse/strumpack/master> relative tolerance for compression 37229e0a805SPieter Ghysels 37329e0a805SPieter Ghysels Logically Collective 37429e0a805SPieter Ghysels 37529e0a805SPieter Ghysels Input Parameters: 37629e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` 37729e0a805SPieter Ghysels 37829e0a805SPieter Ghysels Output Parameter: 37929e0a805SPieter Ghysels . rtol - relative compression tolerance 38029e0a805SPieter Ghysels 38129e0a805SPieter Ghysels Level: intermediate 38229e0a805SPieter Ghysels 3831d27aa22SBarry Smith .seealso: `MATSOLVERSTRUMPACK`, [](ch_matrices), `Mat`, `MatGetFactor()`, `MatSTRUMPACKSetCompRelTol()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()` 38429e0a805SPieter Ghysels @*/ 38529e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetCompRelTol(Mat F, PetscReal *rtol) 38629e0a805SPieter Ghysels { 38729e0a805SPieter Ghysels PetscFunctionBegin; 38829e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 38929e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKGetCompRelTol_C", (Mat, PetscReal *), (F, rtol)); 39029e0a805SPieter Ghysels PetscValidLogicalCollectiveReal(F, *rtol, 2); 3913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 39234f43fa5SPieter Ghysels } 39334f43fa5SPieter Ghysels 39429e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetCompAbsTol_STRUMPACK(Mat F, PetscReal atol) 395d71ae5a4SJacob Faibussowitsch { 39635e8bcc9SJunchao Zhang STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 39734f43fa5SPieter Ghysels 39834f43fa5SPieter Ghysels PetscFunctionBegin; 39929e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_compression_abs_tol", STRUMPACK_set_compression_abs_tol(*S, atol)); 40029e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 40129e0a805SPieter Ghysels } 40229e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetCompAbsTol_STRUMPACK(Mat F, PetscReal *atol) 40329e0a805SPieter Ghysels { 40429e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 40529e0a805SPieter Ghysels 40629e0a805SPieter Ghysels PetscFunctionBegin; 40729e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression_abs_tol", *atol = (PetscReal)STRUMPACK_compression_abs_tol(*S)); 4083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 40934f43fa5SPieter Ghysels } 410e5a36eccSStefano Zampini 41134f43fa5SPieter Ghysels /*@ 4121d27aa22SBarry Smith MatSTRUMPACKSetCompAbsTol - Set STRUMPACK <https://portal.nersc.gov/project/sparse/strumpack/master> absolute tolerance for compression 413147403d9SBarry Smith 414c3339decSBarry Smith Logically Collective 41534f43fa5SPieter Ghysels 41634f43fa5SPieter Ghysels Input Parameters: 4172ef1f0ffSBarry Smith + F - the factored matrix obtained by calling `MatGetFactor()` 41834f43fa5SPieter Ghysels - atol - absolute compression tolerance 41934f43fa5SPieter Ghysels 42011a5261eSBarry Smith Options Database Key: 42129e0a805SPieter Ghysels . -mat_strumpack_compression_abs_tol <1e-10> - Absolute compression tolerance, when using `-pctype ilu` 42234f43fa5SPieter Ghysels 42329e0a805SPieter Ghysels Level: intermediate 42434f43fa5SPieter Ghysels 4251d27aa22SBarry Smith .seealso: `MATSOLVERSTRUMPACK`, [](ch_matrices), `Mat`, `MatGetFactor()`, `MatSTRUMPACKGetCompAbsTol()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()` 42634f43fa5SPieter Ghysels @*/ 42729e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetCompAbsTol(Mat F, PetscReal atol) 428d71ae5a4SJacob Faibussowitsch { 42934f43fa5SPieter Ghysels PetscFunctionBegin; 43034f43fa5SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 43134f43fa5SPieter Ghysels PetscValidLogicalCollectiveReal(F, atol, 2); 43229e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKSetCompAbsTol_C", (Mat, PetscReal), (F, atol)); 4333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 43434f43fa5SPieter Ghysels } 43534f43fa5SPieter Ghysels /*@ 4361d27aa22SBarry Smith MatSTRUMPACKGetCompAbsTol - Get STRUMPACK <https://portal.nersc.gov/project/sparse/strumpack/master> absolute tolerance for compression 437147403d9SBarry Smith 438c3339decSBarry Smith Logically Collective 43934f43fa5SPieter Ghysels 44034f43fa5SPieter Ghysels Input Parameters: 44129e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` 44234f43fa5SPieter Ghysels 44329e0a805SPieter Ghysels Output Parameter: 44429e0a805SPieter Ghysels . atol - absolute compression tolerance 44534f43fa5SPieter Ghysels 44629e0a805SPieter Ghysels Level: intermediate 44734f43fa5SPieter Ghysels 4481d27aa22SBarry Smith .seealso: `MATSOLVERSTRUMPACK`, [](ch_matrices), `Mat`, `MatGetFactor()`, `MatSTRUMPACKSetCompAbsTol()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()` 44934f43fa5SPieter Ghysels @*/ 45029e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetCompAbsTol(Mat F, PetscReal *atol) 451d71ae5a4SJacob Faibussowitsch { 45234f43fa5SPieter Ghysels PetscFunctionBegin; 45334f43fa5SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 45429e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKGetCompAbsTol_C", (Mat, PetscReal *), (F, atol)); 45529e0a805SPieter Ghysels PetscValidLogicalCollectiveReal(F, *atol, 2); 4563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 45734f43fa5SPieter Ghysels } 45834f43fa5SPieter Ghysels 45929e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetCompLeafSize_STRUMPACK(Mat F, PetscInt leaf_size) 460d71ae5a4SJacob Faibussowitsch { 46135e8bcc9SJunchao Zhang STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 462a36bf211SPieter Ghysels 463a36bf211SPieter Ghysels PetscFunctionBegin; 46429e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_compression_leaf_size", STRUMPACK_set_compression_leaf_size(*S, leaf_size)); 46529e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 46629e0a805SPieter Ghysels } 46729e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetCompLeafSize_STRUMPACK(Mat F, PetscInt *leaf_size) 46829e0a805SPieter Ghysels { 46929e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 47029e0a805SPieter Ghysels 47129e0a805SPieter Ghysels PetscFunctionBegin; 47229e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression_leaf_size", *leaf_size = (PetscInt)STRUMPACK_compression_leaf_size(*S)); 4733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 474a36bf211SPieter Ghysels } 475e5a36eccSStefano Zampini 476a36bf211SPieter Ghysels /*@ 4771d27aa22SBarry Smith MatSTRUMPACKSetCompLeafSize - Set STRUMPACK <https://portal.nersc.gov/project/sparse/strumpack/master> leaf size for HSS, BLR, HODLR... 478147403d9SBarry Smith 479c3339decSBarry Smith Logically Collective 480a36bf211SPieter Ghysels 481a36bf211SPieter Ghysels Input Parameters: 48229e0a805SPieter Ghysels + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 48329e0a805SPieter Ghysels - leaf_size - Size of diagonal blocks in rank-structured approximation 484a36bf211SPieter Ghysels 48511a5261eSBarry Smith Options Database Key: 48629e0a805SPieter Ghysels . -mat_strumpack_compression_leaf_size - Size of diagonal blocks in rank-structured approximation, when using `-pctype ilu` 487a36bf211SPieter Ghysels 48829e0a805SPieter Ghysels Level: intermediate 489a36bf211SPieter Ghysels 4901d27aa22SBarry Smith .seealso: `MATSOLVERSTRUMPACK`, [](ch_matrices), `Mat`, `MatGetFactor()`, `MatSTRUMPACKGetCompLeafSize()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()` 491a36bf211SPieter Ghysels @*/ 49229e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetCompLeafSize(Mat F, PetscInt leaf_size) 493d71ae5a4SJacob Faibussowitsch { 494a36bf211SPieter Ghysels PetscFunctionBegin; 495a36bf211SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 496a36bf211SPieter Ghysels PetscValidLogicalCollectiveInt(F, leaf_size, 2); 49729e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKSetCompLeafSize_C", (Mat, PetscInt), (F, leaf_size)); 4983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 499a36bf211SPieter Ghysels } 500291d0ed5SPieter Ghysels /*@ 5011d27aa22SBarry Smith MatSTRUMPACKGetCompLeafSize - Get STRUMPACK <https://portal.nersc.gov/project/sparse/strumpack/master> leaf size for HSS, BLR, HODLR... 502147403d9SBarry Smith 503c3339decSBarry Smith Logically Collective 504291d0ed5SPieter Ghysels 505291d0ed5SPieter Ghysels Input Parameters: 50629e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 507291d0ed5SPieter Ghysels 50829e0a805SPieter Ghysels Output Parameter: 50929e0a805SPieter Ghysels . leaf_size - Size of diagonal blocks in rank-structured approximation 510291d0ed5SPieter Ghysels 51129e0a805SPieter Ghysels Level: intermediate 512291d0ed5SPieter Ghysels 5131d27aa22SBarry Smith .seealso: `MATSOLVERSTRUMPACK`, [](ch_matrices), `Mat`, `MatGetFactor()`, `MatSTRUMPACKSetCompLeafSize()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()` 514291d0ed5SPieter Ghysels @*/ 51529e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetCompLeafSize(Mat F, PetscInt *leaf_size) 516d71ae5a4SJacob Faibussowitsch { 517291d0ed5SPieter Ghysels PetscFunctionBegin; 518291d0ed5SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 51929e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKGetCompLeafSize_C", (Mat, PetscInt *), (F, leaf_size)); 52029e0a805SPieter Ghysels PetscValidLogicalCollectiveInt(F, *leaf_size, 2); 52129e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 52229e0a805SPieter Ghysels } 52329e0a805SPieter Ghysels 52429e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetGeometricNxyz_STRUMPACK(Mat F, PetscInt nx, PetscInt ny, PetscInt nz) 52529e0a805SPieter Ghysels { 52629e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 52729e0a805SPieter Ghysels 52829e0a805SPieter Ghysels PetscFunctionBegin; 52929e0a805SPieter Ghysels if (nx < 1) { 53029e0a805SPieter Ghysels PetscCheck(nx == PETSC_DECIDE || nx == PETSC_DEFAULT, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_OUTOFRANGE, "nx < 1"); 53129e0a805SPieter Ghysels nx = 1; 53229e0a805SPieter Ghysels } 53329e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_nx", STRUMPACK_set_nx(*S, nx)); 53429e0a805SPieter Ghysels if (ny < 1) { 53529e0a805SPieter Ghysels PetscCheck(ny == PETSC_DECIDE || ny == PETSC_DEFAULT, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_OUTOFRANGE, "ny < 1"); 53629e0a805SPieter Ghysels ny = 1; 53729e0a805SPieter Ghysels } 53829e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_ny", STRUMPACK_set_ny(*S, ny)); 53929e0a805SPieter Ghysels if (nz < 1) { 54029e0a805SPieter Ghysels PetscCheck(nz == PETSC_DECIDE || nz == PETSC_DEFAULT, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_OUTOFRANGE, "nz < 1"); 54129e0a805SPieter Ghysels nz = 1; 54229e0a805SPieter Ghysels } 54329e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_nz", STRUMPACK_set_nz(*S, nz)); 54429e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 54529e0a805SPieter Ghysels } 54629e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetGeometricComponents_STRUMPACK(Mat F, PetscInt nc) 54729e0a805SPieter Ghysels { 54829e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 54929e0a805SPieter Ghysels 55029e0a805SPieter Ghysels PetscFunctionBegin; 55129e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_components", STRUMPACK_set_components(*S, nc)); 55229e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 55329e0a805SPieter Ghysels } 55429e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetGeometricWidth_STRUMPACK(Mat F, PetscInt w) 55529e0a805SPieter Ghysels { 55629e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 55729e0a805SPieter Ghysels 55829e0a805SPieter Ghysels PetscFunctionBegin; 55929e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_separator_width", STRUMPACK_set_separator_width(*S, w)); 56029e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 56129e0a805SPieter Ghysels } 56229e0a805SPieter Ghysels 56329e0a805SPieter Ghysels /*@ 5641d27aa22SBarry Smith MatSTRUMPACKSetGeometricNxyz - Set STRUMPACK <https://portal.nersc.gov/project/sparse/strumpack/master> mesh x, y and z dimensions, for use with GEOMETRIC ordering. 56529e0a805SPieter Ghysels 56629e0a805SPieter Ghysels Logically Collective 56729e0a805SPieter Ghysels 56829e0a805SPieter Ghysels Input Parameters: 56929e0a805SPieter Ghysels + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 57029e0a805SPieter Ghysels . nx - x dimension of the mesh 57129e0a805SPieter Ghysels . ny - y dimension of the mesh 57229e0a805SPieter Ghysels - nz - z dimension of the mesh 57329e0a805SPieter Ghysels 57429e0a805SPieter Ghysels Level: intermediate 57529e0a805SPieter Ghysels 5761d27aa22SBarry Smith Note: 57729e0a805SPieter Ghysels If the mesh is two (or one) dimensional one can use 1, `PETSC_DECIDE` or `PETSC_DEFAULT` 57829e0a805SPieter Ghysels for the missing z (and y) dimensions. 57929e0a805SPieter Ghysels 5801d27aa22SBarry Smith .seealso: `MATSOLVERSTRUMPACK`, `MatGetFactor()` 58129e0a805SPieter Ghysels @*/ 58229e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetGeometricNxyz(Mat F, PetscInt nx, PetscInt ny, PetscInt nz) 58329e0a805SPieter Ghysels { 58429e0a805SPieter Ghysels PetscFunctionBegin; 58529e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 58629e0a805SPieter Ghysels PetscValidLogicalCollectiveInt(F, nx, 2); 58729e0a805SPieter Ghysels PetscValidLogicalCollectiveInt(F, ny, 3); 58829e0a805SPieter Ghysels PetscValidLogicalCollectiveInt(F, nz, 4); 58929e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKSetGeometricNxyz_C", (Mat, PetscInt, PetscInt, PetscInt), (F, nx, ny, nz)); 59029e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 59129e0a805SPieter Ghysels } 59229e0a805SPieter Ghysels /*@ 5931d27aa22SBarry Smith MatSTRUMPACKSetGeometricComponents - Set STRUMPACK <https://portal.nersc.gov/project/sparse/strumpack/master> 5941d27aa22SBarry Smith number of degrees of freedom per mesh point, for use with GEOMETRIC ordering. 59529e0a805SPieter Ghysels 59629e0a805SPieter Ghysels Logically Collective 59729e0a805SPieter Ghysels 59829e0a805SPieter Ghysels Input Parameters: 59929e0a805SPieter Ghysels + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 60029e0a805SPieter Ghysels - nc - Number of components/dof's per grid point 60129e0a805SPieter Ghysels 60229e0a805SPieter Ghysels Options Database Key: 60329e0a805SPieter Ghysels . -mat_strumpack_geometric_components <1> - Number of components per mesh point, for geometric nested dissection ordering 60429e0a805SPieter Ghysels 60529e0a805SPieter Ghysels Level: intermediate 60629e0a805SPieter Ghysels 6071d27aa22SBarry Smith .seealso: `MATSOLVERSTRUMPACK`, `MatGetFactor()` 60829e0a805SPieter Ghysels @*/ 60929e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetGeometricComponents(Mat F, PetscInt nc) 61029e0a805SPieter Ghysels { 61129e0a805SPieter Ghysels PetscFunctionBegin; 61229e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 61329e0a805SPieter Ghysels PetscValidLogicalCollectiveInt(F, nc, 2); 61429e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKSetGeometricComponents_C", (Mat, PetscInt), (F, nc)); 61529e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 61629e0a805SPieter Ghysels } 61729e0a805SPieter Ghysels /*@ 6181d27aa22SBarry Smith MatSTRUMPACKSetGeometricWidth - Set STRUMPACK <https://portal.nersc.gov/project/sparse/strumpack/master> width of the separator, for use with GEOMETRIC ordering. 61929e0a805SPieter Ghysels 62029e0a805SPieter Ghysels Logically Collective 62129e0a805SPieter Ghysels 62229e0a805SPieter Ghysels Input Parameters: 62329e0a805SPieter Ghysels + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 62429e0a805SPieter Ghysels - w - width of the separator 62529e0a805SPieter Ghysels 62629e0a805SPieter Ghysels Options Database Key: 62729e0a805SPieter Ghysels . -mat_strumpack_geometric_width <1> - Width of the separator of the mesh, for geometric nested dissection ordering 62829e0a805SPieter Ghysels 62929e0a805SPieter Ghysels Level: intermediate 63029e0a805SPieter Ghysels 6311d27aa22SBarry Smith .seealso: `MATSOLVERSTRUMPACK`, `MatGetFactor()` 63229e0a805SPieter Ghysels @*/ 63329e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetGeometricWidth(Mat F, PetscInt w) 63429e0a805SPieter Ghysels { 63529e0a805SPieter Ghysels PetscFunctionBegin; 63629e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 63729e0a805SPieter Ghysels PetscValidLogicalCollectiveInt(F, w, 2); 63829e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKSetGeometricWidth_C", (Mat, PetscInt), (F, w)); 63929e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 64029e0a805SPieter Ghysels } 64129e0a805SPieter Ghysels 64229e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetCompMinSepSize_STRUMPACK(Mat F, PetscInt min_sep_size) 64329e0a805SPieter Ghysels { 64429e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 64529e0a805SPieter Ghysels 64629e0a805SPieter Ghysels PetscFunctionBegin; 64729e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_compression_min_sep_size", STRUMPACK_set_compression_min_sep_size(*S, min_sep_size)); 64829e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 64929e0a805SPieter Ghysels } 65029e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetCompMinSepSize_STRUMPACK(Mat F, PetscInt *min_sep_size) 65129e0a805SPieter Ghysels { 65229e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 65329e0a805SPieter Ghysels 65429e0a805SPieter Ghysels PetscFunctionBegin; 65529e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression_min_sep_size", *min_sep_size = (PetscInt)STRUMPACK_compression_min_sep_size(*S)); 65629e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 65729e0a805SPieter Ghysels } 65829e0a805SPieter Ghysels 65929e0a805SPieter Ghysels /*@ 6601d27aa22SBarry Smith MatSTRUMPACKSetCompMinSepSize - Set STRUMPACK <https://portal.nersc.gov/project/sparse/strumpack/master> minimum separator size for low-rank approximation 66129e0a805SPieter Ghysels 66229e0a805SPieter Ghysels Logically Collective 66329e0a805SPieter Ghysels 66429e0a805SPieter Ghysels Input Parameters: 66529e0a805SPieter Ghysels + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 66629e0a805SPieter Ghysels - min_sep_size - minimum dense matrix size for low-rank approximation 66729e0a805SPieter Ghysels 66829e0a805SPieter Ghysels Options Database Key: 66929e0a805SPieter Ghysels . -mat_strumpack_compression_min_sep_size <min_sep_size> - Minimum size of dense sub-block for low-rank compression 67029e0a805SPieter Ghysels 67129e0a805SPieter Ghysels Level: intermediate 67229e0a805SPieter Ghysels 6731d27aa22SBarry Smith .seealso: `MATSOLVERSTRUMPACK`, `MatGetFactor()`, `MatSTRUMPACKGetCompMinSepSize()` 67429e0a805SPieter Ghysels @*/ 67529e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetCompMinSepSize(Mat F, PetscInt min_sep_size) 67629e0a805SPieter Ghysels { 67729e0a805SPieter Ghysels PetscFunctionBegin; 67829e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 67929e0a805SPieter Ghysels PetscValidLogicalCollectiveInt(F, min_sep_size, 2); 68029e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKSetCompMinSepSize_C", (Mat, PetscInt), (F, min_sep_size)); 68129e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 68229e0a805SPieter Ghysels } 68329e0a805SPieter Ghysels /*@ 6841d27aa22SBarry Smith MatSTRUMPACKGetCompMinSepSize - Get STRUMPACK <https://portal.nersc.gov/project/sparse/strumpack/master> minimum separator size for low-rank approximation 68529e0a805SPieter Ghysels 68629e0a805SPieter Ghysels Logically Collective 68729e0a805SPieter Ghysels 68829e0a805SPieter Ghysels Input Parameters: 68929e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 69029e0a805SPieter Ghysels 69129e0a805SPieter Ghysels Output Parameter: 69229e0a805SPieter Ghysels . min_sep_size - minimum dense matrix size for low-rank approximation 69329e0a805SPieter Ghysels 69429e0a805SPieter Ghysels Level: intermediate 69529e0a805SPieter Ghysels 6961d27aa22SBarry Smith .seealso: `MATSOLVERSTRUMPACK`, `MatGetFactor()`, `MatSTRUMPACKSetCompMinSepSize()` 69729e0a805SPieter Ghysels @*/ 69829e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetCompMinSepSize(Mat F, PetscInt *min_sep_size) 69929e0a805SPieter Ghysels { 70029e0a805SPieter Ghysels PetscFunctionBegin; 70129e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 70229e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKGetCompMinSepSize_C", (Mat, PetscInt *), (F, min_sep_size)); 70329e0a805SPieter Ghysels PetscValidLogicalCollectiveInt(F, *min_sep_size, 2); 70429e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 70529e0a805SPieter Ghysels } 70629e0a805SPieter Ghysels 70729e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetCompLossyPrecision_STRUMPACK(Mat F, PetscInt lossy_prec) 70829e0a805SPieter Ghysels { 70929e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 71029e0a805SPieter Ghysels 71129e0a805SPieter Ghysels PetscFunctionBegin; 71229e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_compression_lossy_precision", STRUMPACK_set_compression_lossy_precision(*S, lossy_prec)); 71329e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 71429e0a805SPieter Ghysels } 71529e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetCompLossyPrecision_STRUMPACK(Mat F, PetscInt *lossy_prec) 71629e0a805SPieter Ghysels { 71729e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 71829e0a805SPieter Ghysels 71929e0a805SPieter Ghysels PetscFunctionBegin; 72029e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression_lossy_precision", *lossy_prec = (PetscInt)STRUMPACK_compression_lossy_precision(*S)); 72129e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 72229e0a805SPieter Ghysels } 72329e0a805SPieter Ghysels 72429e0a805SPieter Ghysels /*@ 7251d27aa22SBarry Smith MatSTRUMPACKSetCompLossyPrecision - Set STRUMPACK <https://portal.nersc.gov/project/sparse/strumpack/master> precision for lossy compression (requires ZFP support) 72629e0a805SPieter Ghysels 72729e0a805SPieter Ghysels Logically Collective 72829e0a805SPieter Ghysels 72929e0a805SPieter Ghysels Input Parameters: 73029e0a805SPieter Ghysels + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 73129e0a805SPieter Ghysels - lossy_prec - Number of bitplanes to use in lossy compression 73229e0a805SPieter Ghysels 73329e0a805SPieter Ghysels Options Database Key: 73429e0a805SPieter Ghysels . -mat_strumpack_compression_lossy_precision <lossy_prec> - Precision when using lossy compression [1-64], when using `-pctype ilu -mat_strumpack_compression MAT_STRUMPACK_COMPRESSION_TYPE_LOSSY` 73529e0a805SPieter Ghysels 73629e0a805SPieter Ghysels Level: intermediate 73729e0a805SPieter Ghysels 7381d27aa22SBarry Smith .seealso: `MATSOLVERSTRUMPACK`, `MatGetFactor()`, `MatSTRUMPACKGetCompLossyPrecision()` 73929e0a805SPieter Ghysels @*/ 74029e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetCompLossyPrecision(Mat F, PetscInt lossy_prec) 74129e0a805SPieter Ghysels { 74229e0a805SPieter Ghysels PetscFunctionBegin; 74329e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 74429e0a805SPieter Ghysels PetscValidLogicalCollectiveInt(F, lossy_prec, 2); 74529e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKSetCompLossyPrecision_C", (Mat, PetscInt), (F, lossy_prec)); 74629e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 74729e0a805SPieter Ghysels } 74829e0a805SPieter Ghysels /*@ 7491d27aa22SBarry Smith MatSTRUMPACKGetCompLossyPrecision - Get STRUMPACK <https://portal.nersc.gov/project/sparse/strumpack/master> precision for lossy compression (requires ZFP support) 75029e0a805SPieter Ghysels 75129e0a805SPieter Ghysels Logically Collective 75229e0a805SPieter Ghysels 75329e0a805SPieter Ghysels Input Parameters: 75429e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 75529e0a805SPieter Ghysels 75629e0a805SPieter Ghysels Output Parameter: 75729e0a805SPieter Ghysels . lossy_prec - Number of bitplanes to use in lossy compression 75829e0a805SPieter Ghysels 75929e0a805SPieter Ghysels Level: intermediate 76029e0a805SPieter Ghysels 7611d27aa22SBarry Smith .seealso: `MATSOLVERSTRUMPACK`, `MatGetFactor()`, `MatSTRUMPACKSetCompLossyPrecision()` 76229e0a805SPieter Ghysels @*/ 76329e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetCompLossyPrecision(Mat F, PetscInt *lossy_prec) 76429e0a805SPieter Ghysels { 76529e0a805SPieter Ghysels PetscFunctionBegin; 76629e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 76729e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKGetCompLossyPrecision_C", (Mat, PetscInt *), (F, lossy_prec)); 76829e0a805SPieter Ghysels PetscValidLogicalCollectiveInt(F, *lossy_prec, 2); 76929e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 77029e0a805SPieter Ghysels } 77129e0a805SPieter Ghysels 77229e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetCompButterflyLevels_STRUMPACK(Mat F, PetscInt bfly_lvls) 77329e0a805SPieter Ghysels { 77429e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 77529e0a805SPieter Ghysels 77629e0a805SPieter Ghysels PetscFunctionBegin; 77729e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_compression_butterfly_levels", STRUMPACK_set_compression_butterfly_levels(*S, bfly_lvls)); 77829e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 77929e0a805SPieter Ghysels } 78029e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetCompButterflyLevels_STRUMPACK(Mat F, PetscInt *bfly_lvls) 78129e0a805SPieter Ghysels { 78229e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 78329e0a805SPieter Ghysels 78429e0a805SPieter Ghysels PetscFunctionBegin; 78529e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression_butterfly_levels", *bfly_lvls = (PetscInt)STRUMPACK_compression_butterfly_levels(*S)); 78629e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 78729e0a805SPieter Ghysels } 78829e0a805SPieter Ghysels 78929e0a805SPieter Ghysels /*@ 7901d27aa22SBarry Smith MatSTRUMPACKSetCompButterflyLevels - Set STRUMPACK <https://portal.nersc.gov/project/sparse/strumpack/master> 7911d27aa22SBarry Smith number of butterfly levels in HODLR compression (requires ButterflyPACK support) 79229e0a805SPieter Ghysels 79329e0a805SPieter Ghysels Logically Collective 79429e0a805SPieter Ghysels 79529e0a805SPieter Ghysels Input Parameters: 79629e0a805SPieter Ghysels + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 79729e0a805SPieter Ghysels - bfly_lvls - Number of levels of butterfly compression in HODLR compression 79829e0a805SPieter Ghysels 79929e0a805SPieter Ghysels Options Database Key: 8001d27aa22SBarry Smith . -mat_strumpack_compression_butterfly_levels <bfly_lvls> - Number of levels in the hierarchically off-diagonal matrix for which to use butterfly, 8011d27aa22SBarry Smith when using `-pctype ilu`, (BLR_)HODLR compression 80229e0a805SPieter Ghysels 80329e0a805SPieter Ghysels Level: intermediate 80429e0a805SPieter Ghysels 8051d27aa22SBarry Smith .seealso: `MATSOLVERSTRUMPACK`, `MatGetFactor()`, `MatSTRUMPACKGetCompButterflyLevels()` 80629e0a805SPieter Ghysels @*/ 80729e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetCompButterflyLevels(Mat F, PetscInt bfly_lvls) 80829e0a805SPieter Ghysels { 80929e0a805SPieter Ghysels PetscFunctionBegin; 81029e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 81129e0a805SPieter Ghysels PetscValidLogicalCollectiveInt(F, bfly_lvls, 2); 81229e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKSetButterflyLevels_C", (Mat, PetscInt), (F, bfly_lvls)); 81329e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 81429e0a805SPieter Ghysels } 81529e0a805SPieter Ghysels /*@ 8161d27aa22SBarry Smith MatSTRUMPACKGetCompButterflyLevels - Get STRUMPACK <https://portal.nersc.gov/project/sparse/strumpack/master> 8171d27aa22SBarry Smith number of butterfly levels in HODLR compression (requires ButterflyPACK support) 81829e0a805SPieter Ghysels 81929e0a805SPieter Ghysels Logically Collective 82029e0a805SPieter Ghysels 82129e0a805SPieter Ghysels Input Parameters: 82229e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 82329e0a805SPieter Ghysels 82429e0a805SPieter Ghysels Output Parameter: 82529e0a805SPieter Ghysels . bfly_lvls - Number of levels of butterfly compression in HODLR compression 82629e0a805SPieter Ghysels 82729e0a805SPieter Ghysels Level: intermediate 82829e0a805SPieter Ghysels 8291d27aa22SBarry Smith .seealso: `MATSOLVERSTRUMPACK`, `MatGetFactor()`, `MatSTRUMPACKSetCompButterflyLevels()` 83029e0a805SPieter Ghysels @*/ 83129e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetCompButterflyLevels(Mat F, PetscInt *bfly_lvls) 83229e0a805SPieter Ghysels { 83329e0a805SPieter Ghysels PetscFunctionBegin; 83429e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 83529e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKGetButterflyLevels_C", (Mat, PetscInt *), (F, bfly_lvls)); 83629e0a805SPieter Ghysels PetscValidLogicalCollectiveInt(F, *bfly_lvls, 2); 8373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 838291d0ed5SPieter Ghysels } 839291d0ed5SPieter Ghysels 840d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolve_STRUMPACK(Mat A, Vec b_mpi, Vec x) 841d71ae5a4SJacob Faibussowitsch { 84235e8bcc9SJunchao Zhang STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)A->data; 8437d6ea485SPieter Ghysels STRUMPACK_RETURN_CODE sp_err; 8447d6ea485SPieter Ghysels const PetscScalar *bptr; 8457d6ea485SPieter Ghysels PetscScalar *xptr; 8467d6ea485SPieter Ghysels 8477d6ea485SPieter Ghysels PetscFunctionBegin; 8489566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &xptr)); 8499566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(b_mpi, &bptr)); 8500d08a34dSPieter Ghysels 851e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_solve", sp_err = STRUMPACK_solve(*S, (PetscScalar *)bptr, xptr, 0)); 8520d08a34dSPieter Ghysels switch (sp_err) { 853d71ae5a4SJacob Faibussowitsch case STRUMPACK_SUCCESS: 854d71ae5a4SJacob Faibussowitsch break; 8559371c9d4SSatish Balay case STRUMPACK_MATRIX_NOT_SET: { 8569371c9d4SSatish Balay SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix was not set"); 8579371c9d4SSatish Balay break; 8589371c9d4SSatish Balay } 8599371c9d4SSatish Balay case STRUMPACK_REORDERING_ERROR: { 8609371c9d4SSatish Balay SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix reordering failed"); 8619371c9d4SSatish Balay break; 8629371c9d4SSatish Balay } 863d71ae5a4SJacob Faibussowitsch default: 864d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: solve failed"); 8657d6ea485SPieter Ghysels } 8669566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &xptr)); 8679566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(b_mpi, &bptr)); 8683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8697d6ea485SPieter Ghysels } 8707d6ea485SPieter Ghysels 871d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMatSolve_STRUMPACK(Mat A, Mat B_mpi, Mat X) 872d71ae5a4SJacob Faibussowitsch { 87329e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)A->data; 87429e0a805SPieter Ghysels STRUMPACK_RETURN_CODE sp_err; 8757d6ea485SPieter Ghysels PetscBool flg; 87629e0a805SPieter Ghysels PetscInt m = A->rmap->n, nrhs; 87729e0a805SPieter Ghysels const PetscScalar *bptr; 87829e0a805SPieter Ghysels PetscScalar *xptr; 8797d6ea485SPieter Ghysels 8807d6ea485SPieter Ghysels PetscFunctionBegin; 8819566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)B_mpi, &flg, MATSEQDENSE, MATMPIDENSE, NULL)); 8825f80ce2aSJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Matrix B must be MATDENSE matrix"); 8839566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)X, &flg, MATSEQDENSE, MATMPIDENSE, NULL)); 8845f80ce2aSJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Matrix X must be MATDENSE matrix"); 88529e0a805SPieter Ghysels 88629e0a805SPieter Ghysels PetscCall(MatGetSize(B_mpi, NULL, &nrhs)); 88729e0a805SPieter Ghysels PetscCall(MatDenseGetArray(X, &xptr)); 88829e0a805SPieter Ghysels PetscCall(MatDenseGetArrayRead(B_mpi, &bptr)); 88929e0a805SPieter Ghysels 89029e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_matsolve", sp_err = STRUMPACK_matsolve(*S, nrhs, bptr, m, xptr, m, 0)); 89129e0a805SPieter Ghysels switch (sp_err) { 89229e0a805SPieter Ghysels case STRUMPACK_SUCCESS: 89329e0a805SPieter Ghysels break; 89429e0a805SPieter Ghysels case STRUMPACK_MATRIX_NOT_SET: { 89529e0a805SPieter Ghysels SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix was not set"); 89629e0a805SPieter Ghysels break; 89729e0a805SPieter Ghysels } 89829e0a805SPieter Ghysels case STRUMPACK_REORDERING_ERROR: { 89929e0a805SPieter Ghysels SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix reordering failed"); 90029e0a805SPieter Ghysels break; 90129e0a805SPieter Ghysels } 90229e0a805SPieter Ghysels default: 90329e0a805SPieter Ghysels SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: solve failed"); 90429e0a805SPieter Ghysels } 90529e0a805SPieter Ghysels PetscCall(MatDenseRestoreArrayRead(B_mpi, &bptr)); 90629e0a805SPieter Ghysels PetscCall(MatDenseRestoreArray(X, &xptr)); 9073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9087d6ea485SPieter Ghysels } 9097d6ea485SPieter Ghysels 910d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_Info_STRUMPACK(Mat A, PetscViewer viewer) 911d71ae5a4SJacob Faibussowitsch { 912ad0c5e61SPieter Ghysels PetscFunctionBegin; 913ad0c5e61SPieter Ghysels /* check if matrix is strumpack type */ 9143ba16761SJacob Faibussowitsch if (A->ops->solve != MatSolve_STRUMPACK) PetscFunctionReturn(PETSC_SUCCESS); 9159566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "STRUMPACK sparse solver!\n")); 9163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 917ad0c5e61SPieter Ghysels } 918ad0c5e61SPieter Ghysels 919d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_STRUMPACK(Mat A, PetscViewer viewer) 920d71ae5a4SJacob Faibussowitsch { 9219f196a02SMartin Diehl PetscBool isascii; 922ad0c5e61SPieter Ghysels PetscViewerFormat format; 923ad0c5e61SPieter Ghysels 924ad0c5e61SPieter Ghysels PetscFunctionBegin; 9259f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 9269f196a02SMartin Diehl if (isascii) { 9279566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 9289566063dSJacob Faibussowitsch if (format == PETSC_VIEWER_ASCII_INFO) PetscCall(MatView_Info_STRUMPACK(A, viewer)); 929ad0c5e61SPieter Ghysels } 9303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 931ad0c5e61SPieter Ghysels } 9327d6ea485SPieter Ghysels 933d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLUFactorNumeric_STRUMPACK(Mat F, Mat A, const MatFactorInfo *info) 934d71ae5a4SJacob Faibussowitsch { 93535e8bcc9SJunchao Zhang STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 9367d6ea485SPieter Ghysels STRUMPACK_RETURN_CODE sp_err; 93729e0a805SPieter Ghysels Mat Aloc; 93829e0a805SPieter Ghysels const PetscScalar *av; 93929e0a805SPieter Ghysels const PetscInt *ai = NULL, *aj = NULL; 94029e0a805SPieter Ghysels PetscInt M = A->rmap->N, m = A->rmap->n, dummy; 94129e0a805SPieter Ghysels PetscBool ismpiaij, isseqaij, flg; 9427d6ea485SPieter Ghysels 9437d6ea485SPieter Ghysels PetscFunctionBegin; 94429e0a805SPieter Ghysels PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATSEQAIJ, &isseqaij)); 94529e0a805SPieter Ghysels PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATMPIAIJ, &ismpiaij)); 94629e0a805SPieter Ghysels if (ismpiaij) { 94729e0a805SPieter Ghysels PetscCall(MatMPIAIJGetLocalMat(A, MAT_INITIAL_MATRIX, &Aloc)); 94829e0a805SPieter Ghysels } else if (isseqaij) { 94929e0a805SPieter Ghysels PetscCall(PetscObjectReference((PetscObject)A)); 95029e0a805SPieter Ghysels Aloc = A; 95129e0a805SPieter Ghysels } else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not for type %s", ((PetscObject)A)->type_name); 95229e0a805SPieter Ghysels 95329e0a805SPieter Ghysels PetscCall(MatGetRowIJ(Aloc, 0, PETSC_FALSE, PETSC_FALSE, &dummy, &ai, &aj, &flg)); 95429e0a805SPieter Ghysels PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "GetRowIJ failed"); 95529e0a805SPieter Ghysels PetscCall(MatSeqAIJGetArrayRead(Aloc, &av)); 95629e0a805SPieter Ghysels 95729e0a805SPieter Ghysels if (ismpiaij) { 95829e0a805SPieter Ghysels const PetscInt *dist = NULL; 95929e0a805SPieter Ghysels PetscCall(MatGetOwnershipRanges(A, &dist)); 96029e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_distributed_csr_matrix", STRUMPACK_set_distributed_csr_matrix(*S, &m, ai, aj, av, dist, 0)); 96129e0a805SPieter Ghysels } else if (isseqaij) { 96229e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_csr_matrix", STRUMPACK_set_csr_matrix(*S, &M, ai, aj, av, 0)); 9637d6ea485SPieter Ghysels } 9647d6ea485SPieter Ghysels 96529e0a805SPieter Ghysels PetscCall(MatRestoreRowIJ(Aloc, 0, PETSC_FALSE, PETSC_FALSE, &dummy, &ai, &aj, &flg)); 96629e0a805SPieter Ghysels PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "RestoreRowIJ failed"); 96729e0a805SPieter Ghysels PetscCall(MatSeqAIJRestoreArrayRead(Aloc, &av)); 96829e0a805SPieter Ghysels PetscCall(MatDestroy(&Aloc)); 96929e0a805SPieter Ghysels 9707d6ea485SPieter Ghysels /* Reorder and Factor the matrix. */ 9717d6ea485SPieter Ghysels /* TODO figure out how to avoid reorder if the matrix values changed, but the pattern remains the same. */ 972e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_reorder", sp_err = STRUMPACK_reorder(*S)); 973e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_factor", sp_err = STRUMPACK_factor(*S)); 9740d08a34dSPieter Ghysels switch (sp_err) { 975d71ae5a4SJacob Faibussowitsch case STRUMPACK_SUCCESS: 976d71ae5a4SJacob Faibussowitsch break; 9779371c9d4SSatish Balay case STRUMPACK_MATRIX_NOT_SET: { 9789371c9d4SSatish Balay SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix was not set"); 9799371c9d4SSatish Balay break; 9809371c9d4SSatish Balay } 9819371c9d4SSatish Balay case STRUMPACK_REORDERING_ERROR: { 9829371c9d4SSatish Balay SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix reordering failed"); 9839371c9d4SSatish Balay break; 9849371c9d4SSatish Balay } 985d71ae5a4SJacob Faibussowitsch default: 986d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: factorization failed"); 9877d6ea485SPieter Ghysels } 988cb250fa3SPieter Ghysels F->assembled = PETSC_TRUE; 989cb250fa3SPieter Ghysels F->preallocated = PETSC_TRUE; 9903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9917d6ea485SPieter Ghysels } 9927d6ea485SPieter Ghysels 993d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLUFactorSymbolic_STRUMPACK(Mat F, Mat A, IS r, IS c, const MatFactorInfo *info) 994d71ae5a4SJacob Faibussowitsch { 9957d6ea485SPieter Ghysels PetscFunctionBegin; 9967d6ea485SPieter Ghysels F->ops->lufactornumeric = MatLUFactorNumeric_STRUMPACK; 9977d6ea485SPieter Ghysels F->ops->solve = MatSolve_STRUMPACK; 9987d6ea485SPieter Ghysels F->ops->matsolve = MatMatSolve_STRUMPACK; 9993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10007d6ea485SPieter Ghysels } 10017d6ea485SPieter Ghysels 1002d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorGetSolverType_aij_strumpack(Mat A, MatSolverType *type) 1003d71ae5a4SJacob Faibussowitsch { 10047d6ea485SPieter Ghysels PetscFunctionBegin; 10057d6ea485SPieter Ghysels *type = MATSOLVERSTRUMPACK; 10063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10077d6ea485SPieter Ghysels } 10087d6ea485SPieter Ghysels 1009575704cbSPieter Ghysels /*MC 101029e0a805SPieter Ghysels MATSOLVERSTRUMPACK = "strumpack" - A solver package providing a direct sparse solver (PCLU) 10111d27aa22SBarry Smith and a preconditioner (PCILU) using low-rank compression via the external package STRUMPACK <https://portal.nersc.gov/project/sparse/strumpack/master>. 1012575704cbSPieter Ghysels 101329e0a805SPieter Ghysels Use `./configure --download-strumpack --download-metis` to have PETSc installed with STRUMPACK. 101429e0a805SPieter Ghysels 101529e0a805SPieter Ghysels For full functionality, add `--download-slate --download-magma --download-parmetis --download-ptscotch --download-zfp --download-butterflypack`. 101629e0a805SPieter Ghysels SLATE provides GPU support in the multi-GPU setting, providing ScaLAPACK functionality but with GPU acceleration. 101729e0a805SPieter Ghysels MAGMA can optionally be used for on node GPU support instead cuBLAS/cuSOLVER, and performs slightly better. 101829e0a805SPieter Ghysels ParMETIS and PTScotch can be used for parallel fill-reducing ordering. 101929e0a805SPieter Ghysels ZFP is used for floating point compression of the sparse factors (LOSSY or LOSSLESS compression). 102029e0a805SPieter Ghysels ButterflyPACK is used for HODLR (Hierarchically Off-Diagonal Low Rank) and HODBF (Hierarchically Off-Diagonal Butterfly) compression of the sparse factors. 102129e0a805SPieter Ghysels 1022575704cbSPieter Ghysels Options Database Keys: 102329e0a805SPieter Ghysels + -mat_strumpack_verbose - Enable verbose output 102429e0a805SPieter Ghysels . -mat_strumpack_compression - Type of rank-structured compression in sparse LU factors (choose one of) NONE HSS BLR HODLR BLR_HODLR ZFP_BLR_HODLR LOSSLESS LOSSY 102529e0a805SPieter Ghysels . -mat_strumpack_compression_rel_tol - Relative compression tolerance, when using `-pctype ilu` 102629e0a805SPieter Ghysels . -mat_strumpack_compression_abs_tol - Absolute compression tolerance, when using `-pctype ilu` 102729e0a805SPieter Ghysels . -mat_strumpack_compression_min_sep_size - Minimum size of separator for rank-structured compression, when using `-pctype ilu` 102829e0a805SPieter Ghysels . -mat_strumpack_compression_leaf_size - Size of diagonal blocks in rank-structured approximation, when using `-pctype ilu` 102929e0a805SPieter Ghysels . -mat_strumpack_compression_lossy_precision - Precision when using lossy compression [1-64], when using `-pctype ilu`, compression LOSSY (requires ZFP support) 103029e0a805SPieter Ghysels . -mat_strumpack_compression_butterfly_levels - Number of levels in the hierarchically off-diagonal matrix for which to use butterfly, when using `-pctype ilu`, (BLR_)HODLR compression (requires ButterflyPACK support) 103129e0a805SPieter Ghysels . -mat_strumpack_gpu - Enable GPU acceleration in numerical factorization (not supported for all compression types) 103229e0a805SPieter Ghysels . -mat_strumpack_colperm <TRUE> - Permute matrix to make diagonal nonzeros 103329e0a805SPieter Ghysels . -mat_strumpack_reordering <METIS> - Sparsity reducing matrix reordering (choose one of) NATURAL METIS PARMETIS SCOTCH PTSCOTCH RCM GEOMETRIC AMD MMD AND MLF SPECTRAL 103429e0a805SPieter Ghysels . -mat_strumpack_geometric_xyz <1,1,1> - Mesh x,y,z dimensions, for use with GEOMETRIC ordering 103529e0a805SPieter Ghysels . -mat_strumpack_geometric_components <1> - Number of components per mesh point, for geometric nested dissection ordering 103629e0a805SPieter Ghysels . -mat_strumpack_geometric_width <1> - Width of the separator of the mesh, for geometric nested dissection ordering 103729e0a805SPieter Ghysels - -mat_strumpack_metis_nodeNDP - Use METIS_NodeNDP instead of METIS_NodeND, for a more balanced tree 1038575704cbSPieter Ghysels 1039575704cbSPieter Ghysels Level: beginner 1040575704cbSPieter Ghysels 10411d27aa22SBarry Smith Notes: 10421d27aa22SBarry Smith Recommended use is 1 MPI process per GPU. 10431d27aa22SBarry Smith 10441d27aa22SBarry Smith Use `-pc_type lu` `-pc_factor_mat_solver_type strumpack` to use this as an exact (direct) solver. 10451d27aa22SBarry Smith 10461d27aa22SBarry Smith Use `-pc_type ilu` `-pc_factor_mat_solver_type strumpack` to enable low-rank compression (i.e, use as a preconditioner), by default using block low rank (BLR). 10471d27aa22SBarry Smith 10481d27aa22SBarry Smith Works with `MATAIJ` matrices 10491d27aa22SBarry Smith 105029e0a805SPieter Ghysels HODLR, BLR_HODBF and ZFP_BLR_HODLR compression require STRUMPACK to be configured with ButterflyPACK support (`--download-butterflypack`). 105129e0a805SPieter Ghysels 105229e0a805SPieter Ghysels LOSSY, LOSSLESS and ZFP_BLR_HODLR compression require STRUMPACK to be configured with ZFP support (`--download-zfp`). 105329e0a805SPieter Ghysels 10541d27aa22SBarry Smith .seealso: `MATSOLVERSTRUMPACK`, [](ch_matrices), `Mat`, `PCLU`, `PCILU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`, `PCFactorSetMatSolverType()`, `MatSolverType`, 105529e0a805SPieter Ghysels `MatGetFactor()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKReordering`, `MatSTRUMPACKCompressionType`, `MatSTRUMPACKSetColPerm()`. 1056575704cbSPieter Ghysels M*/ 1057d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetFactor_aij_strumpack(Mat A, MatFactorType ftype, Mat *F) 1058d71ae5a4SJacob Faibussowitsch { 10597d6ea485SPieter Ghysels Mat B; 10607d6ea485SPieter Ghysels PetscInt M = A->rmap->N, N = A->cmap->N; 106129e0a805SPieter Ghysels PetscBool verb, flg, set; 106229e0a805SPieter Ghysels PetscReal ctol; 106329e0a805SPieter Ghysels PetscInt min_sep_size, leaf_size, nxyz[3], nrdims, nc, w; 106429e0a805SPieter Ghysels #if defined(STRUMPACK_USE_ZFP) 106529e0a805SPieter Ghysels PetscInt lossy_prec; 106629e0a805SPieter Ghysels #endif 106729e0a805SPieter Ghysels #if defined(STRUMPACK_USE_BPACK) 106829e0a805SPieter Ghysels PetscInt bfly_lvls; 106929e0a805SPieter Ghysels #endif 107029e0a805SPieter Ghysels #if defined(STRUMPACK_USE_SLATE_SCALAPACK) 107129e0a805SPieter Ghysels PetscMPIInt mpithreads; 107229e0a805SPieter Ghysels #endif 107334f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S; 107434f43fa5SPieter Ghysels STRUMPACK_INTERFACE iface; 107529e0a805SPieter Ghysels STRUMPACK_REORDERING_STRATEGY ndcurrent, ndvalue; 107629e0a805SPieter Ghysels STRUMPACK_COMPRESSION_TYPE compcurrent, compvalue; 10779371c9d4SSatish Balay const STRUMPACK_PRECISION table[2][2][2] = { 10789371c9d4SSatish Balay {{STRUMPACK_FLOATCOMPLEX_64, STRUMPACK_DOUBLECOMPLEX_64}, {STRUMPACK_FLOAT_64, STRUMPACK_DOUBLE_64}}, 10799371c9d4SSatish Balay {{STRUMPACK_FLOATCOMPLEX, STRUMPACK_DOUBLECOMPLEX}, {STRUMPACK_FLOAT, STRUMPACK_DOUBLE} } 10809371c9d4SSatish Balay }; 10814ac6704cSBarry Smith const STRUMPACK_PRECISION prec = table[(sizeof(PetscInt) == 8) ? 0 : 1][(PETSC_SCALAR == PETSC_COMPLEX) ? 0 : 1][(PETSC_REAL == PETSC_FLOAT) ? 0 : 1]; 108229e0a805SPieter Ghysels const char *const STRUMPACKNDTypes[] = {"NATURAL", "METIS", "PARMETIS", "SCOTCH", "PTSCOTCH", "RCM", "GEOMETRIC", "AMD", "MMD", "AND", "MLF", "SPECTRAL", "STRUMPACKNDTypes", "", 0}; 108329e0a805SPieter Ghysels const char *const CompTypes[] = {"NONE", "HSS", "BLR", "HODLR", "BLR_HODLR", "ZFP_BLR_HODLR", "LOSSLESS", "LOSSY", "CompTypes", "", 0}; 10847d6ea485SPieter Ghysels 10857d6ea485SPieter Ghysels PetscFunctionBegin; 108629e0a805SPieter Ghysels #if defined(STRUMPACK_USE_SLATE_SCALAPACK) 108729e0a805SPieter Ghysels PetscCallMPI(MPI_Query_thread(&mpithreads)); 108829e0a805SPieter Ghysels PetscCheck(mpithreads == MPI_THREAD_MULTIPLE, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP_SYS, "SLATE requires MPI_THREAD_MULTIPLE"); 108929e0a805SPieter Ghysels #endif 10907d6ea485SPieter Ghysels /* Create the factorization matrix */ 10919566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B)); 10929566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, M, N)); 109335e8bcc9SJunchao Zhang PetscCall(PetscStrallocpy("strumpack", &((PetscObject)B)->type_name)); 109435e8bcc9SJunchao Zhang PetscCall(MatSetUp(B)); 109529e0a805SPieter Ghysels PetscCall(MatSeqAIJSetPreallocation(B, 0, NULL)); 109629e0a805SPieter Ghysels PetscCall(MatMPIAIJSetPreallocation(B, 0, NULL, 0, NULL)); 109766e17bc3SBarry Smith B->trivialsymbolic = PETSC_TRUE; 1098*966bd95aSPierre Jolivet PetscCheck(ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU, PETSC_COMM_SELF, PETSC_ERR_SUP, "Factor type not supported"); 10997d6ea485SPieter Ghysels B->ops->lufactorsymbolic = MatLUFactorSymbolic_STRUMPACK; 1100575704cbSPieter Ghysels B->ops->ilufactorsymbolic = MatLUFactorSymbolic_STRUMPACK; 110135e8bcc9SJunchao Zhang B->ops->getinfo = MatGetInfo_External; 11027d6ea485SPieter Ghysels B->ops->view = MatView_STRUMPACK; 11037d6ea485SPieter Ghysels B->ops->destroy = MatDestroy_STRUMPACK; 11047d6ea485SPieter Ghysels B->ops->getdiagonal = MatGetDiagonal_STRUMPACK; 11059566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_aij_strumpack)); 11069566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetReordering_C", MatSTRUMPACKSetReordering_STRUMPACK)); 110729e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetReordering_C", MatSTRUMPACKGetReordering_STRUMPACK)); 11089566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetColPerm_C", MatSTRUMPACKSetColPerm_STRUMPACK)); 110929e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetColPerm_C", MatSTRUMPACKGetColPerm_STRUMPACK)); 111029e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetGeometricNxyz_C", MatSTRUMPACKSetGeometricNxyz_STRUMPACK)); 111129e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetGeometricComponents_C", MatSTRUMPACKSetGeometricComponents_STRUMPACK)); 111229e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetGeometricWidth_C", MatSTRUMPACKSetGeometricWidth_STRUMPACK)); 111329e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetGPU_C", MatSTRUMPACKSetGPU_STRUMPACK)); 111429e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetGPU_C", MatSTRUMPACKGetGPU_STRUMPACK)); 111529e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetCompression_C", MatSTRUMPACKSetCompression_STRUMPACK)); 111629e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetCompression_C", MatSTRUMPACKGetCompression_STRUMPACK)); 111729e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetCompRelTol_C", MatSTRUMPACKSetCompRelTol_STRUMPACK)); 111829e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetCompRelTol_C", MatSTRUMPACKGetCompRelTol_STRUMPACK)); 111929e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetCompAbsTol_C", MatSTRUMPACKSetCompAbsTol_STRUMPACK)); 112029e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetCompAbsTol_C", MatSTRUMPACKGetCompAbsTol_STRUMPACK)); 112129e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetCompLeafSize_C", MatSTRUMPACKSetCompLeafSize_STRUMPACK)); 112229e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetCompLeafSize_C", MatSTRUMPACKGetCompLeafSize_STRUMPACK)); 112329e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetCompMinSepSize_C", MatSTRUMPACKSetCompMinSepSize_STRUMPACK)); 112429e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetCompMinSepSize_C", MatSTRUMPACKGetCompMinSepSize_STRUMPACK)); 112529e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetCompLossyPrecision_C", MatSTRUMPACKSetCompLossyPrecision_STRUMPACK)); 112629e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetCompLossyPrecision_C", MatSTRUMPACKGetCompLossyPrecision_STRUMPACK)); 112729e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetCompButterflyLevels_C", MatSTRUMPACKSetCompButterflyLevels_STRUMPACK)); 112829e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetCompButterflyLevels_C", MatSTRUMPACKGetCompButterflyLevels_STRUMPACK)); 1129575704cbSPieter Ghysels B->factortype = ftype; 113029e0a805SPieter Ghysels 113129e0a805SPieter Ghysels /* set solvertype */ 113229e0a805SPieter Ghysels PetscCall(PetscFree(B->solvertype)); 113329e0a805SPieter Ghysels PetscCall(PetscStrallocpy(MATSOLVERSTRUMPACK, &B->solvertype)); 113429e0a805SPieter Ghysels 11354dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&S)); 113635e8bcc9SJunchao Zhang B->data = S; 11370d08a34dSPieter Ghysels 113835e8bcc9SJunchao Zhang PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATSEQAIJ, &flg)); /* A might be MATSEQAIJCUSPARSE */ 11390d08a34dSPieter Ghysels iface = flg ? STRUMPACK_MT : STRUMPACK_MPI_DIST; 11407d6ea485SPieter Ghysels 114126cc229bSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)B), ((PetscObject)B)->prefix, "STRUMPACK Options", "Mat"); 1142575704cbSPieter Ghysels verb = PetscLogPrintInfo ? PETSC_TRUE : PETSC_FALSE; 11439566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_strumpack_verbose", "Print STRUMPACK information", "None", verb, &verb, NULL)); 11447d6ea485SPieter Ghysels 1145e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_init", STRUMPACK_init(S, PetscObjectComm((PetscObject)A), prec, iface, 0, NULL, verb)); 114655c022e5SPieter Ghysels 114729e0a805SPieter Ghysels /* By default, no compression is done. Compression is enabled when the user enables it with */ 114829e0a805SPieter Ghysels /* -mat_strumpack_compression with anything else than NONE, or when selecting ilu */ 114929e0a805SPieter Ghysels /* preconditioning, in which case we default to STRUMPACK_BLR compression. */ 115029e0a805SPieter Ghysels /* When compression is enabled, the STRUMPACK solver becomes an incomplete */ 115188382b05SPieter Ghysels /* (or approximate) LU factorization. */ 115229e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression", compcurrent = STRUMPACK_compression(*S)); 115329e0a805SPieter Ghysels PetscCall(PetscOptionsEnum("-mat_strumpack_compression", "Rank-structured compression type", "None", CompTypes, (PetscEnum)compcurrent, (PetscEnum *)&compvalue, &set)); 115429e0a805SPieter Ghysels if (set) { 115529e0a805SPieter Ghysels PetscCall(MatSTRUMPACKSetCompression(B, (MatSTRUMPACKCompressionType)compvalue)); 115629e0a805SPieter Ghysels } else { 115729e0a805SPieter Ghysels if (ftype == MAT_FACTOR_ILU) { PetscStackCallExternalVoid("STRUMPACK_set_compression", STRUMPACK_set_compression(*S, STRUMPACK_BLR)); } 115888382b05SPieter Ghysels } 115955c022e5SPieter Ghysels 116029e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression_rel_tol", ctol = (PetscReal)STRUMPACK_compression_rel_tol(*S)); 116129e0a805SPieter Ghysels PetscCall(PetscOptionsReal("-mat_strumpack_compression_rel_tol", "Relative compression tolerance", "None", ctol, &ctol, &set)); 116229e0a805SPieter Ghysels if (set) PetscStackCallExternalVoid("STRUMPACK_set_compression_rel_tol", STRUMPACK_set_compression_rel_tol(*S, (double)ctol)); 116329e0a805SPieter Ghysels 116429e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression_abs_tol", ctol = (PetscReal)STRUMPACK_compression_abs_tol(*S)); 116529e0a805SPieter Ghysels PetscCall(PetscOptionsReal("-mat_strumpack_compression_abs_tol", "Absolute compression tolerance", "None", ctol, &ctol, &set)); 116629e0a805SPieter Ghysels if (set) PetscStackCallExternalVoid("STRUMPACK_set_compression_abs_tol", STRUMPACK_set_compression_abs_tol(*S, (double)ctol)); 116729e0a805SPieter Ghysels 116829e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression_min_sep_size", min_sep_size = (PetscInt)STRUMPACK_compression_min_sep_size(*S)); 116929e0a805SPieter Ghysels PetscCall(PetscOptionsInt("-mat_strumpack_compression_min_sep_size", "Minimum size of separator for compression", "None", min_sep_size, &min_sep_size, &set)); 117029e0a805SPieter Ghysels if (set) PetscStackCallExternalVoid("STRUMPACK_set_compression_min_sep_size", STRUMPACK_set_compression_min_sep_size(*S, (int)min_sep_size)); 117129e0a805SPieter Ghysels 117229e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression_leaf_size", leaf_size = (PetscInt)STRUMPACK_compression_leaf_size(*S)); 117329e0a805SPieter Ghysels PetscCall(PetscOptionsInt("-mat_strumpack_compression_leaf_size", "Size of diagonal blocks in rank-structured approximation", "None", leaf_size, &leaf_size, &set)); 117429e0a805SPieter Ghysels if (set) PetscStackCallExternalVoid("STRUMPACK_set_compression_leaf_size", STRUMPACK_set_compression_leaf_size(*S, (int)leaf_size)); 117529e0a805SPieter Ghysels 117629e0a805SPieter Ghysels #if defined(STRUMPACK_USE_ZFP) 117729e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression_lossy_precision", lossy_prec = (PetscInt)STRUMPACK_compression_lossy_precision(*S)); 117829e0a805SPieter Ghysels PetscCall(PetscOptionsInt("-mat_strumpack_compression_lossy_precision", "Number of bitplanes to use in lossy compression", "None", lossy_prec, &lossy_prec, &set)); 117929e0a805SPieter Ghysels if (set) PetscStackCallExternalVoid("STRUMPACK_set_compression_lossy_precision", STRUMPACK_set_compression_lossy_precision(*S, (int)lossy_prec)); 118029e0a805SPieter Ghysels #endif 118129e0a805SPieter Ghysels 118229e0a805SPieter Ghysels #if defined(STRUMPACK_USE_BPACK) 118329e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression_butterfly_levels", bfly_lvls = (PetscInt)STRUMPACK_compression_butterfly_levels(*S)); 118429e0a805SPieter Ghysels PetscCall(PetscOptionsInt("-mat_strumpack_compression_butterfly_levels", "Number of levels in the HODLR matrix for which to use butterfly compression", "None", bfly_lvls, &bfly_lvls, &set)); 118529e0a805SPieter Ghysels if (set) PetscStackCallExternalVoid("STRUMPACK_set_compression_butterfly_levels", STRUMPACK_set_compression_butterfly_levels(*S, (int)bfly_lvls)); 118629e0a805SPieter Ghysels #endif 118729e0a805SPieter Ghysels 118829e0a805SPieter Ghysels #if defined(STRUMPACK_USE_CUDA) || defined(STRUMPACK_USE_HIP) || defined(STRUMPACK_USE_SYCL) 118929e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_use_gpu", flg = (STRUMPACK_use_gpu(*S) == 0) ? PETSC_FALSE : PETSC_TRUE); 119029e0a805SPieter Ghysels PetscCall(PetscOptionsBool("-mat_strumpack_gpu", "Enable GPU acceleration (not supported for all compression types)", "None", flg, &flg, &set)); 119129e0a805SPieter Ghysels if (set) MatSTRUMPACKSetGPU(B, flg); 119229e0a805SPieter Ghysels #endif 119329e0a805SPieter Ghysels 119429e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_matching", flg = (STRUMPACK_matching(*S) == STRUMPACK_MATCHING_NONE) ? PETSC_FALSE : PETSC_TRUE); 119529e0a805SPieter Ghysels PetscCall(PetscOptionsBool("-mat_strumpack_colperm", "Find a col perm to get nonzero diagonal", "None", flg, &flg, &set)); 119629e0a805SPieter Ghysels if (set) PetscStackCallExternalVoid("STRUMPACK_set_matching", STRUMPACK_set_matching(*S, flg ? STRUMPACK_MATCHING_MAX_DIAGONAL_PRODUCT_SCALING : STRUMPACK_MATCHING_NONE)); 119729e0a805SPieter Ghysels 119829e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_reordering_method", ndcurrent = STRUMPACK_reordering_method(*S)); 119929e0a805SPieter Ghysels PetscCall(PetscOptionsEnum("-mat_strumpack_reordering", "Sparsity reducing matrix reordering", "None", STRUMPACKNDTypes, (PetscEnum)ndcurrent, (PetscEnum *)&ndvalue, &set)); 120029e0a805SPieter Ghysels if (set) PetscStackCallExternalVoid("STRUMPACK_set_reordering_method", STRUMPACK_set_reordering_method(*S, ndvalue)); 120129e0a805SPieter Ghysels 120229e0a805SPieter Ghysels /* geometric ordering, for a regular 1D/2D/3D mesh in the natural ordering, */ 120329e0a805SPieter Ghysels /* with nc DOF's per gridpoint, and possibly a wider stencil */ 120429e0a805SPieter Ghysels nrdims = 3; 120529e0a805SPieter Ghysels nxyz[0] = nxyz[1] = nxyz[2] = 1; 120629e0a805SPieter Ghysels PetscCall(PetscOptionsIntArray("-mat_strumpack_geometric_xyz", "Mesh sizes nx,ny,nz (Use 1 for default)", "", nxyz, &nrdims, &set)); 120729e0a805SPieter Ghysels if (set) { 120829e0a805SPieter Ghysels PetscCheck(nrdims == 1 || nrdims == 2 || nrdims == 3, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_OUTOFRANGE, "'-mat_strumpack_geometric_xyz' requires 1, 2, or 3 values."); 120929e0a805SPieter Ghysels PetscCall(MatSTRUMPACKSetGeometricNxyz(B, (int)nxyz[0], (int)nxyz[1], (int)nxyz[2])); 121029e0a805SPieter Ghysels } 121129e0a805SPieter Ghysels PetscCall(PetscOptionsInt("-mat_strumpack_geometric_components", "Number of components per mesh point, for geometric nested dissection ordering", "None", 1, &nc, &set)); 121229e0a805SPieter Ghysels if (set) PetscStackCallExternalVoid("STRUMPACK_set_components", STRUMPACK_set_components(*S, (int)nc)); 121329e0a805SPieter Ghysels PetscCall(PetscOptionsInt("-mat_strumpack_geometric_width", "Width of the separator (for instance a 1D 3-point wide stencil needs a 1 point wide separator, a 1D 5-point stencil needs a 2 point wide separator), for geometric nested dissection ordering", "None", 1, &w, &set)); 121429e0a805SPieter Ghysels if (set) PetscStackCallExternalVoid("STRUMPACK_set_separator_width", STRUMPACK_set_separator_width(*S, (int)w)); 121529e0a805SPieter Ghysels 121629e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_use_METIS_NodeNDP", flg = (STRUMPACK_use_METIS_NodeNDP(*S) == 0) ? PETSC_FALSE : PETSC_TRUE); 121729e0a805SPieter Ghysels PetscCall(PetscOptionsBool("-mat_strumpack_metis_nodeNDP", "Use METIS_NodeNDP instead of METIS_NodeND, for a more balanced tree", "None", flg, &flg, &set)); 121829e0a805SPieter Ghysels if (set) { 121929e0a805SPieter Ghysels if (flg) { 122029e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_enable_METIS_NodeNDP", STRUMPACK_enable_METIS_NodeNDP(*S)); 122129e0a805SPieter Ghysels } else { 122229e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_disable_METIS_NodeNDP", STRUMPACK_disable_METIS_NodeNDP(*S)); 122329e0a805SPieter Ghysels } 122429e0a805SPieter Ghysels } 122529e0a805SPieter Ghysels 122629e0a805SPieter Ghysels /* Disable the outer iterative solver from STRUMPACK. */ 122729e0a805SPieter Ghysels /* When STRUMPACK is used as a direct solver, it will by default do iterative refinement. */ 122829e0a805SPieter Ghysels /* When STRUMPACK is used as an approximate factorization preconditioner (by enabling */ 122929e0a805SPieter Ghysels /* low-rank compression), it will use it's own preconditioned GMRES. Here we can disable */ 123029e0a805SPieter Ghysels /* the outer iterative solver, as PETSc uses STRUMPACK from within a KSP. */ 123129e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_Krylov_solver", STRUMPACK_set_Krylov_solver(*S, STRUMPACK_DIRECT)); 123229e0a805SPieter Ghysels 123329e0a805SPieter Ghysels PetscOptionsEnd(); 123426cc229bSBarry Smith 12357d6ea485SPieter Ghysels *F = B; 12363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12377d6ea485SPieter Ghysels } 12387d6ea485SPieter Ghysels 1239d1f0640dSPierre Jolivet PETSC_INTERN PetscErrorCode MatSolverTypeRegister_STRUMPACK(void) 1240d71ae5a4SJacob Faibussowitsch { 12417d6ea485SPieter Ghysels PetscFunctionBegin; 12429566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSTRUMPACK, MATMPIAIJ, MAT_FACTOR_LU, MatGetFactor_aij_strumpack)); 12439566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSTRUMPACK, MATSEQAIJ, MAT_FACTOR_LU, MatGetFactor_aij_strumpack)); 12449566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSTRUMPACK, MATMPIAIJ, MAT_FACTOR_ILU, MatGetFactor_aij_strumpack)); 12459566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSTRUMPACK, MATSEQAIJ, MAT_FACTOR_ILU, MatGetFactor_aij_strumpack)); 12463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12477d6ea485SPieter Ghysels } 1248