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)); 48575704cbSPieter Ghysels 493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 50575704cbSPieter Ghysels } 51575704cbSPieter Ghysels 52d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSTRUMPACKSetReordering_STRUMPACK(Mat F, MatSTRUMPACKReordering reordering) 53d71ae5a4SJacob Faibussowitsch { 5435e8bcc9SJunchao Zhang STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 55575704cbSPieter Ghysels 56575704cbSPieter Ghysels PetscFunctionBegin; 5729e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_reordering_method", STRUMPACK_set_reordering_method(*S, (STRUMPACK_REORDERING_STRATEGY)reordering)); 5829e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 5929e0a805SPieter Ghysels } 6029e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetReordering_STRUMPACK(Mat F, MatSTRUMPACKReordering *reordering) 6129e0a805SPieter Ghysels { 6229e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 6329e0a805SPieter Ghysels 6429e0a805SPieter Ghysels PetscFunctionBegin; 6529e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_reordering_method", *reordering = (MatSTRUMPACKReordering)STRUMPACK_reordering_method(*S)); 663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6734f43fa5SPieter Ghysels } 68e5a36eccSStefano Zampini 69542aee0fSPieter Ghysels /*@ 7034f43fa5SPieter Ghysels MatSTRUMPACKSetReordering - Set STRUMPACK fill-reducing reordering 7134f43fa5SPieter Ghysels 7229e0a805SPieter Ghysels Logically Collective 7329e0a805SPieter Ghysels 7434f43fa5SPieter Ghysels Input Parameters: 7529e0a805SPieter Ghysels + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 7629e0a805SPieter Ghysels - reordering - the code to be used to find the fill-reducing reordering 7734f43fa5SPieter Ghysels 7811a5261eSBarry Smith Options Database Key: 792ef1f0ffSBarry Smith . -mat_strumpack_reordering <METIS> - Sparsity reducing matrix reordering, see `MatSTRUMPACKReordering` 8034f43fa5SPieter Ghysels 8129e0a805SPieter Ghysels Level: intermediate 8234f43fa5SPieter Ghysels 8334f43fa5SPieter Ghysels References: 8429e0a805SPieter Ghysels . * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/ 8534f43fa5SPieter Ghysels 8629e0a805SPieter Ghysels .seealso: [](ch_matrices), `Mat`, `MatSTRUMPACKReordering`, `MatGetFactor()`, `MatSTRUMPACKSetColPerm()`, `MatSTRUMPACKGetReordering()` 87542aee0fSPieter Ghysels @*/ 88d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSTRUMPACKSetReordering(Mat F, MatSTRUMPACKReordering reordering) 89d71ae5a4SJacob Faibussowitsch { 9034f43fa5SPieter Ghysels PetscFunctionBegin; 9134f43fa5SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 9234f43fa5SPieter Ghysels PetscValidLogicalCollectiveEnum(F, reordering, 2); 93cac4c232SBarry Smith PetscTryMethod(F, "MatSTRUMPACKSetReordering_C", (Mat, MatSTRUMPACKReordering), (F, reordering)); 943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 95575704cbSPieter Ghysels } 9629e0a805SPieter Ghysels /*@ 9729e0a805SPieter Ghysels MatSTRUMPACKGetReordering - Get STRUMPACK fill-reducing reordering 9829e0a805SPieter Ghysels 9929e0a805SPieter Ghysels Logically Collective 10029e0a805SPieter Ghysels 10129e0a805SPieter Ghysels Input Parameters: 10229e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 10329e0a805SPieter Ghysels 10429e0a805SPieter Ghysels Output Parameter: 10529e0a805SPieter Ghysels . reordering - the code to be used to find the fill-reducing reordering 10629e0a805SPieter Ghysels 10729e0a805SPieter Ghysels Level: intermediate 10829e0a805SPieter Ghysels 10929e0a805SPieter Ghysels References: 11029e0a805SPieter Ghysels . * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/ 11129e0a805SPieter Ghysels 11229e0a805SPieter Ghysels .seealso: [](ch_matrices), `Mat`, `MatSTRUMPACKReordering`, `MatGetFactor()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()` 11329e0a805SPieter Ghysels @*/ 11429e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetReordering(Mat F, MatSTRUMPACKReordering *reordering) 11529e0a805SPieter Ghysels { 11629e0a805SPieter Ghysels PetscFunctionBegin; 11729e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 11829e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKGetReordering_C", (Mat, MatSTRUMPACKReordering *), (F, reordering)); 11929e0a805SPieter Ghysels PetscValidLogicalCollectiveEnum(F, *reordering, 2); 12029e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 12129e0a805SPieter Ghysels } 122575704cbSPieter Ghysels 123d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSTRUMPACKSetColPerm_STRUMPACK(Mat F, PetscBool cperm) 124d71ae5a4SJacob Faibussowitsch { 12535e8bcc9SJunchao Zhang STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 12634f43fa5SPieter Ghysels 12734f43fa5SPieter Ghysels PetscFunctionBegin; 12829e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_matching", STRUMPACK_set_matching(*S, cperm ? STRUMPACK_MATCHING_MAX_DIAGONAL_PRODUCT_SCALING : STRUMPACK_MATCHING_NONE)); 12929e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 13029e0a805SPieter Ghysels } 13129e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetColPerm_STRUMPACK(Mat F, PetscBool *cperm) 13229e0a805SPieter Ghysels { 13329e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 13429e0a805SPieter Ghysels 13529e0a805SPieter Ghysels PetscFunctionBegin; 13629e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_matching", *cperm = (PetscBool)(STRUMPACK_matching(*S) != STRUMPACK_MATCHING_NONE)); 1373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13834f43fa5SPieter Ghysels } 139e5a36eccSStefano Zampini 140575704cbSPieter Ghysels /*@ 14134f43fa5SPieter Ghysels MatSTRUMPACKSetColPerm - Set whether STRUMPACK should try to permute the columns of the matrix in order to get a nonzero diagonal 142147403d9SBarry Smith 143c3339decSBarry Smith Logically Collective 144575704cbSPieter Ghysels 145575704cbSPieter Ghysels Input Parameters: 1462ef1f0ffSBarry Smith + F - the factored matrix obtained by calling `MatGetFactor()` 14711a5261eSBarry Smith - cperm - `PETSC_TRUE` to permute (internally) the columns of the matrix 148575704cbSPieter Ghysels 14911a5261eSBarry Smith Options Database Key: 150147403d9SBarry Smith . -mat_strumpack_colperm <cperm> - true to use the permutation 151575704cbSPieter Ghysels 15229e0a805SPieter Ghysels Level: intermediate 153575704cbSPieter Ghysels 154575704cbSPieter Ghysels References: 15529e0a805SPieter Ghysels . * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/ 156575704cbSPieter Ghysels 15729e0a805SPieter Ghysels .seealso: [](ch_matrices), `MatSTRUMPACKSetReordering()`, `Mat`, `MatGetFactor()`, `MatSTRUMPACKGetColPerm()` 158575704cbSPieter Ghysels @*/ 159d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSTRUMPACKSetColPerm(Mat F, PetscBool cperm) 160d71ae5a4SJacob Faibussowitsch { 161575704cbSPieter Ghysels PetscFunctionBegin; 162575704cbSPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 16334f43fa5SPieter Ghysels PetscValidLogicalCollectiveBool(F, cperm, 2); 164cac4c232SBarry Smith PetscTryMethod(F, "MatSTRUMPACKSetColPerm_C", (Mat, PetscBool), (F, cperm)); 1653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 166575704cbSPieter Ghysels } 16729e0a805SPieter Ghysels /*@ 16829e0a805SPieter Ghysels MatSTRUMPACKGetColPerm - Get whether STRUMPACK will try to permute the columns of the matrix in order to get a nonzero diagonal 169575704cbSPieter Ghysels 17029e0a805SPieter Ghysels Logically Collective 17129e0a805SPieter Ghysels 17229e0a805SPieter Ghysels Input Parameters: 17329e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` 17429e0a805SPieter Ghysels 17529e0a805SPieter Ghysels Output Parameter: 17629e0a805SPieter Ghysels . cperm - Indicates whether STRUMPACK will permute columns 17729e0a805SPieter Ghysels 17829e0a805SPieter Ghysels Level: intermediate 17929e0a805SPieter Ghysels 18029e0a805SPieter Ghysels References: 18129e0a805SPieter Ghysels . * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/ 18229e0a805SPieter Ghysels 18329e0a805SPieter Ghysels .seealso: [](ch_matrices), `MatSTRUMPACKSetReordering()`, `Mat`, `MatGetFactor()`, `MatSTRUMPACKSetColPerm()` 18429e0a805SPieter Ghysels @*/ 18529e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetColPerm(Mat F, PetscBool *cperm) 18629e0a805SPieter Ghysels { 18729e0a805SPieter Ghysels PetscFunctionBegin; 18829e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 18929e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKGetColPerm_C", (Mat, PetscBool *), (F, cperm)); 19029e0a805SPieter Ghysels PetscValidLogicalCollectiveBool(F, *cperm, 2); 19129e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 19229e0a805SPieter Ghysels } 19329e0a805SPieter Ghysels 19429e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetGPU_STRUMPACK(Mat F, PetscBool gpu) 195d71ae5a4SJacob Faibussowitsch { 19635e8bcc9SJunchao Zhang STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 19734f43fa5SPieter Ghysels 19834f43fa5SPieter Ghysels PetscFunctionBegin; 19929e0a805SPieter Ghysels if (gpu) { 20029e0a805SPieter Ghysels #if !(defined(STRUMPACK_USE_CUDA) || defined(STRUMPACK_USE_HIP) || defined(STRUMPACK_USE_SYCL)) 20129e0a805SPieter Ghysels PetscCall(PetscPrintf(PETSC_COMM_SELF, " Warning: strumpack was not configured with GPU support\n")); 20229e0a805SPieter Ghysels #endif 20329e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_enable_gpu", STRUMPACK_enable_gpu(*S)); 20429e0a805SPieter Ghysels } else PetscStackCallExternalVoid("STRUMPACK_disable_gpu", STRUMPACK_disable_gpu(*S)); 20529e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 20629e0a805SPieter Ghysels } 20729e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetGPU_STRUMPACK(Mat F, PetscBool *gpu) 20829e0a805SPieter Ghysels { 20929e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 21029e0a805SPieter Ghysels 21129e0a805SPieter Ghysels PetscFunctionBegin; 21229e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_use_gpu", *gpu = (PetscBool)STRUMPACK_use_gpu(*S)); 2133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 21434f43fa5SPieter Ghysels } 215e5a36eccSStefano Zampini 21634f43fa5SPieter Ghysels /*@ 21729e0a805SPieter Ghysels MatSTRUMPACKSetGPU - Set whether STRUMPACK should enable GPU acceleration (not supported for all compression types) 21829e0a805SPieter Ghysels 21929e0a805SPieter Ghysels Logically Collective 22029e0a805SPieter Ghysels 22129e0a805SPieter Ghysels Input Parameters: 22229e0a805SPieter Ghysels + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 22329e0a805SPieter Ghysels - gpu - whether or not to use GPU acceleration 22429e0a805SPieter Ghysels 22529e0a805SPieter Ghysels Options Database Key: 22629e0a805SPieter Ghysels . -mat_strumpack_gpu <gpu> - true to use gpu offload 22729e0a805SPieter Ghysels 22829e0a805SPieter Ghysels Level: intermediate 22929e0a805SPieter Ghysels 23029e0a805SPieter Ghysels References: 23129e0a805SPieter Ghysels . * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/ 23229e0a805SPieter Ghysels 23329e0a805SPieter Ghysels .seealso: `MatGetFactor()`, `MatSTRUMPACKGetGPU()` 23429e0a805SPieter Ghysels @*/ 23529e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetGPU(Mat F, PetscBool gpu) 23629e0a805SPieter Ghysels { 23729e0a805SPieter Ghysels PetscFunctionBegin; 23829e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 23929e0a805SPieter Ghysels PetscValidLogicalCollectiveBool(F, gpu, 2); 24029e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKSetGPU_C", (Mat, PetscBool), (F, gpu)); 24129e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 24229e0a805SPieter Ghysels } 24329e0a805SPieter Ghysels /*@ 24429e0a805SPieter Ghysels MatSTRUMPACKGetGPU - Get whether STRUMPACK will try to use GPU acceleration (not supported for all compression types) 24529e0a805SPieter Ghysels 24629e0a805SPieter Ghysels Logically Collective 24729e0a805SPieter Ghysels 24829e0a805SPieter Ghysels Input Parameters: 24929e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 25029e0a805SPieter Ghysels 25129e0a805SPieter Ghysels Output Parameter: 25229e0a805SPieter Ghysels . gpu - whether or not STRUMPACK will try to use GPU acceleration 25329e0a805SPieter Ghysels 25429e0a805SPieter Ghysels Level: intermediate 25529e0a805SPieter Ghysels 25629e0a805SPieter Ghysels References: 25729e0a805SPieter Ghysels . * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/ 25829e0a805SPieter Ghysels 25929e0a805SPieter Ghysels .seealso: `MatGetFactor()`, `MatSTRUMPACKSetGPU()` 26029e0a805SPieter Ghysels @*/ 26129e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetGPU(Mat F, PetscBool *gpu) 26229e0a805SPieter Ghysels { 26329e0a805SPieter Ghysels PetscFunctionBegin; 26429e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 26529e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKGetGPU_C", (Mat, PetscBool *), (F, gpu)); 26629e0a805SPieter Ghysels PetscValidLogicalCollectiveBool(F, *gpu, 2); 26729e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 26829e0a805SPieter Ghysels } 26929e0a805SPieter Ghysels 27029e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetCompression_STRUMPACK(Mat F, MatSTRUMPACKCompressionType comp) 27129e0a805SPieter Ghysels { 27229e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 27329e0a805SPieter Ghysels 27429e0a805SPieter Ghysels PetscFunctionBegin; 275*c2935fbeSBarry Smith #if !defined(STRUMPACK_USE_BPACK) 27629e0a805SPieter 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"); 27729e0a805SPieter Ghysels #endif 278*c2935fbeSBarry Smith #if !defined(STRUMPACK_USE_ZFP) 27929e0a805SPieter 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"); 28029e0a805SPieter Ghysels #endif 28129e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_compression", STRUMPACK_set_compression(*S, (STRUMPACK_COMPRESSION_TYPE)comp)); 28229e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 28329e0a805SPieter Ghysels } 28429e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetCompression_STRUMPACK(Mat F, MatSTRUMPACKCompressionType *comp) 28529e0a805SPieter Ghysels { 28629e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 28729e0a805SPieter Ghysels 28829e0a805SPieter Ghysels PetscFunctionBegin; 28929e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression", *comp = (MatSTRUMPACKCompressionType)STRUMPACK_compression(*S)); 29029e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 29129e0a805SPieter Ghysels } 29229e0a805SPieter Ghysels 29329e0a805SPieter Ghysels /*@ 29429e0a805SPieter Ghysels MatSTRUMPACKSetCompression - Set STRUMPACK compression type 29529e0a805SPieter Ghysels 29629e0a805SPieter Ghysels Input Parameters: 29729e0a805SPieter Ghysels + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 29829e0a805SPieter Ghysels - comp - Type of compression to be used in the approximate sparse factorization 29929e0a805SPieter Ghysels 30029e0a805SPieter Ghysels Options Database Key: 30129e0a805SPieter 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 30229e0a805SPieter Ghysels 30329e0a805SPieter Ghysels Level: intermediate 30429e0a805SPieter Ghysels 30529e0a805SPieter Ghysels Notes: 30629e0a805SPieter Ghysels Default for `comp` is `MAT_STRUMPACK_COMPRESSION_TYPE_NONE` for `-pc_type lu` and `MAT_STRUMPACK_COMPRESSION_TYPE_BLR` 30729e0a805SPieter Ghysels for `-pc_type ilu` 30829e0a805SPieter Ghysels 30929e0a805SPieter Ghysels References: 31029e0a805SPieter Ghysels . * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/ 31129e0a805SPieter Ghysels 31229e0a805SPieter Ghysels .seealso: `MatGetFactor()`, `MatSTRUMPACKCompressionType`, `MatSTRUMPACKGetCompression()` 31329e0a805SPieter Ghysels @*/ 31429e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetCompression(Mat F, MatSTRUMPACKCompressionType comp) 31529e0a805SPieter Ghysels { 31629e0a805SPieter Ghysels PetscFunctionBegin; 31729e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 31829e0a805SPieter Ghysels PetscValidLogicalCollectiveEnum(F, comp, 2); 31929e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKSetCompression_C", (Mat, MatSTRUMPACKCompressionType), (F, comp)); 32029e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 32129e0a805SPieter Ghysels } 32229e0a805SPieter Ghysels /*@ 32329e0a805SPieter Ghysels MatSTRUMPACKGetCompression - Get STRUMPACK compression type 32429e0a805SPieter Ghysels 32529e0a805SPieter Ghysels Input Parameters: 32629e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 32729e0a805SPieter Ghysels 32829e0a805SPieter Ghysels Output Parameter: 32929e0a805SPieter Ghysels . comp - Type of compression to be used in the approximate sparse factorization 33029e0a805SPieter Ghysels 33129e0a805SPieter Ghysels Level: intermediate 33229e0a805SPieter Ghysels 33329e0a805SPieter Ghysels Notes: 33429e0a805SPieter Ghysels Default is `MAT_STRUMPACK_COMPRESSION_TYPE_NONE` for `-pc_type lu` and `MAT_STRUMPACK_COMPRESSION_TYPE_BLR` for `-pc_type ilu` 33529e0a805SPieter Ghysels 33629e0a805SPieter Ghysels References: 33729e0a805SPieter Ghysels . * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/ 33829e0a805SPieter Ghysels 33929e0a805SPieter Ghysels .seealso: `MatGetFactor()`, `MatSTRUMPACKCompressionType`, `MatSTRUMPACKSetCompression()` 34029e0a805SPieter Ghysels @*/ 34129e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetCompression(Mat F, MatSTRUMPACKCompressionType *comp) 34229e0a805SPieter Ghysels { 34329e0a805SPieter Ghysels PetscFunctionBegin; 34429e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 34529e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKGetCompression_C", (Mat, MatSTRUMPACKCompressionType *), (F, comp)); 34629e0a805SPieter Ghysels PetscValidLogicalCollectiveEnum(F, *comp, 2); 34729e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 34829e0a805SPieter Ghysels } 34929e0a805SPieter Ghysels 35029e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetCompRelTol_STRUMPACK(Mat F, PetscReal rtol) 35129e0a805SPieter Ghysels { 35229e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 35329e0a805SPieter Ghysels 35429e0a805SPieter Ghysels PetscFunctionBegin; 35529e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_compression_rel_tol", STRUMPACK_set_compression_rel_tol(*S, rtol)); 35629e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 35729e0a805SPieter Ghysels } 35829e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetCompRelTol_STRUMPACK(Mat F, PetscReal *rtol) 35929e0a805SPieter Ghysels { 36029e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 36129e0a805SPieter Ghysels 36229e0a805SPieter Ghysels PetscFunctionBegin; 36329e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression_rel_tol", *rtol = (PetscReal)STRUMPACK_compression_rel_tol(*S)); 36429e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 36529e0a805SPieter Ghysels } 36629e0a805SPieter Ghysels 36729e0a805SPieter Ghysels /*@ 36829e0a805SPieter Ghysels MatSTRUMPACKSetCompRelTol - Set STRUMPACK relative tolerance for compression 369147403d9SBarry Smith 370c3339decSBarry Smith Logically Collective 37134f43fa5SPieter Ghysels 37234f43fa5SPieter Ghysels Input Parameters: 3732ef1f0ffSBarry Smith + F - the factored matrix obtained by calling `MatGetFactor()` 37434f43fa5SPieter Ghysels - rtol - relative compression tolerance 37534f43fa5SPieter Ghysels 37611a5261eSBarry Smith Options Database Key: 37729e0a805SPieter Ghysels . -mat_strumpack_compression_rel_tol <1e-4> - Relative compression tolerance, when using `-pctype ilu` 37834f43fa5SPieter Ghysels 37929e0a805SPieter Ghysels Level: intermediate 38034f43fa5SPieter Ghysels 38134f43fa5SPieter Ghysels References: 38229e0a805SPieter Ghysels . * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/ 38334f43fa5SPieter Ghysels 38429e0a805SPieter Ghysels .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatSTRUMPACKGetCompRelTol()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()` 38534f43fa5SPieter Ghysels @*/ 38629e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetCompRelTol(Mat F, PetscReal rtol) 387d71ae5a4SJacob Faibussowitsch { 38834f43fa5SPieter Ghysels PetscFunctionBegin; 38934f43fa5SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 39034f43fa5SPieter Ghysels PetscValidLogicalCollectiveReal(F, rtol, 2); 39129e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKSetCompRelTol_C", (Mat, PetscReal), (F, rtol)); 39229e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 39329e0a805SPieter Ghysels } 39429e0a805SPieter Ghysels /*@ 39529e0a805SPieter Ghysels MatSTRUMPACKGetCompRelTol - Get STRUMPACK relative tolerance for compression 39629e0a805SPieter Ghysels 39729e0a805SPieter Ghysels Logically Collective 39829e0a805SPieter Ghysels 39929e0a805SPieter Ghysels Input Parameters: 40029e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` 40129e0a805SPieter Ghysels 40229e0a805SPieter Ghysels Output Parameter: 40329e0a805SPieter Ghysels . rtol - relative compression tolerance 40429e0a805SPieter Ghysels 40529e0a805SPieter Ghysels Level: intermediate 40629e0a805SPieter Ghysels 40729e0a805SPieter Ghysels References: 40829e0a805SPieter Ghysels . * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/ 40929e0a805SPieter Ghysels 41029e0a805SPieter Ghysels .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatSTRUMPACKSetCompRelTol()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()` 41129e0a805SPieter Ghysels @*/ 41229e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetCompRelTol(Mat F, PetscReal *rtol) 41329e0a805SPieter Ghysels { 41429e0a805SPieter Ghysels PetscFunctionBegin; 41529e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 41629e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKGetCompRelTol_C", (Mat, PetscReal *), (F, rtol)); 41729e0a805SPieter Ghysels PetscValidLogicalCollectiveReal(F, *rtol, 2); 4183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 41934f43fa5SPieter Ghysels } 42034f43fa5SPieter Ghysels 42129e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetCompAbsTol_STRUMPACK(Mat F, PetscReal atol) 422d71ae5a4SJacob Faibussowitsch { 42335e8bcc9SJunchao Zhang STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 42434f43fa5SPieter Ghysels 42534f43fa5SPieter Ghysels PetscFunctionBegin; 42629e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_compression_abs_tol", STRUMPACK_set_compression_abs_tol(*S, atol)); 42729e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 42829e0a805SPieter Ghysels } 42929e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetCompAbsTol_STRUMPACK(Mat F, PetscReal *atol) 43029e0a805SPieter Ghysels { 43129e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 43229e0a805SPieter Ghysels 43329e0a805SPieter Ghysels PetscFunctionBegin; 43429e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression_abs_tol", *atol = (PetscReal)STRUMPACK_compression_abs_tol(*S)); 4353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 43634f43fa5SPieter Ghysels } 437e5a36eccSStefano Zampini 43834f43fa5SPieter Ghysels /*@ 43929e0a805SPieter Ghysels MatSTRUMPACKSetCompAbsTol - Set STRUMPACK absolute tolerance for compression 440147403d9SBarry Smith 441c3339decSBarry Smith Logically Collective 44234f43fa5SPieter Ghysels 44334f43fa5SPieter Ghysels Input Parameters: 4442ef1f0ffSBarry Smith + F - the factored matrix obtained by calling `MatGetFactor()` 44534f43fa5SPieter Ghysels - atol - absolute compression tolerance 44634f43fa5SPieter Ghysels 44711a5261eSBarry Smith Options Database Key: 44829e0a805SPieter Ghysels . -mat_strumpack_compression_abs_tol <1e-10> - Absolute compression tolerance, when using `-pctype ilu` 44934f43fa5SPieter Ghysels 45029e0a805SPieter Ghysels Level: intermediate 45134f43fa5SPieter Ghysels 45234f43fa5SPieter Ghysels References: 45329e0a805SPieter Ghysels . * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/ 45434f43fa5SPieter Ghysels 45529e0a805SPieter Ghysels .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatSTRUMPACKGetCompAbsTol()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()` 45634f43fa5SPieter Ghysels @*/ 45729e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetCompAbsTol(Mat F, PetscReal atol) 458d71ae5a4SJacob Faibussowitsch { 45934f43fa5SPieter Ghysels PetscFunctionBegin; 46034f43fa5SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 46134f43fa5SPieter Ghysels PetscValidLogicalCollectiveReal(F, atol, 2); 46229e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKSetCompAbsTol_C", (Mat, PetscReal), (F, atol)); 4633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 46434f43fa5SPieter Ghysels } 46534f43fa5SPieter Ghysels /*@ 46629e0a805SPieter Ghysels MatSTRUMPACKGetCompAbsTol - Get STRUMPACK absolute tolerance for compression 467147403d9SBarry Smith 468c3339decSBarry Smith Logically Collective 46934f43fa5SPieter Ghysels 47034f43fa5SPieter Ghysels Input Parameters: 47129e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` 47234f43fa5SPieter Ghysels 47329e0a805SPieter Ghysels Output Parameter: 47429e0a805SPieter Ghysels . atol - absolute compression tolerance 47534f43fa5SPieter Ghysels 47629e0a805SPieter Ghysels Level: intermediate 47734f43fa5SPieter Ghysels 47834f43fa5SPieter Ghysels References: 47929e0a805SPieter Ghysels . * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/ 48034f43fa5SPieter Ghysels 48129e0a805SPieter Ghysels .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatSTRUMPACKSetCompAbsTol()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()` 48234f43fa5SPieter Ghysels @*/ 48329e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetCompAbsTol(Mat F, PetscReal *atol) 484d71ae5a4SJacob Faibussowitsch { 48534f43fa5SPieter Ghysels PetscFunctionBegin; 48634f43fa5SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 48729e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKGetCompAbsTol_C", (Mat, PetscReal *), (F, atol)); 48829e0a805SPieter Ghysels PetscValidLogicalCollectiveReal(F, *atol, 2); 4893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 49034f43fa5SPieter Ghysels } 49134f43fa5SPieter Ghysels 49229e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetCompLeafSize_STRUMPACK(Mat F, PetscInt leaf_size) 493d71ae5a4SJacob Faibussowitsch { 49435e8bcc9SJunchao Zhang STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 495a36bf211SPieter Ghysels 496a36bf211SPieter Ghysels PetscFunctionBegin; 49729e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_compression_leaf_size", STRUMPACK_set_compression_leaf_size(*S, leaf_size)); 49829e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 49929e0a805SPieter Ghysels } 50029e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetCompLeafSize_STRUMPACK(Mat F, PetscInt *leaf_size) 50129e0a805SPieter Ghysels { 50229e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 50329e0a805SPieter Ghysels 50429e0a805SPieter Ghysels PetscFunctionBegin; 50529e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression_leaf_size", *leaf_size = (PetscInt)STRUMPACK_compression_leaf_size(*S)); 5063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 507a36bf211SPieter Ghysels } 508e5a36eccSStefano Zampini 509a36bf211SPieter Ghysels /*@ 51029e0a805SPieter Ghysels MatSTRUMPACKSetCompLeafSize - Set STRUMPACK leaf size for HSS, BLR, HODLR... 511147403d9SBarry Smith 512c3339decSBarry Smith Logically Collective 513a36bf211SPieter Ghysels 514a36bf211SPieter Ghysels Input Parameters: 51529e0a805SPieter Ghysels + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 51629e0a805SPieter Ghysels - leaf_size - Size of diagonal blocks in rank-structured approximation 517a36bf211SPieter Ghysels 51811a5261eSBarry Smith Options Database Key: 51929e0a805SPieter Ghysels . -mat_strumpack_compression_leaf_size - Size of diagonal blocks in rank-structured approximation, when using `-pctype ilu` 520a36bf211SPieter Ghysels 52129e0a805SPieter Ghysels Level: intermediate 522a36bf211SPieter Ghysels 523a36bf211SPieter Ghysels References: 52429e0a805SPieter Ghysels . * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/ 525a36bf211SPieter Ghysels 52629e0a805SPieter Ghysels .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatSTRUMPACKGetCompLeafSize()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()` 527a36bf211SPieter Ghysels @*/ 52829e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetCompLeafSize(Mat F, PetscInt leaf_size) 529d71ae5a4SJacob Faibussowitsch { 530a36bf211SPieter Ghysels PetscFunctionBegin; 531a36bf211SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 532a36bf211SPieter Ghysels PetscValidLogicalCollectiveInt(F, leaf_size, 2); 53329e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKSetCompLeafSize_C", (Mat, PetscInt), (F, leaf_size)); 5343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 535a36bf211SPieter Ghysels } 536291d0ed5SPieter Ghysels /*@ 53729e0a805SPieter Ghysels MatSTRUMPACKGetCompLeafSize - Get STRUMPACK leaf size for HSS, BLR, HODLR... 538147403d9SBarry Smith 539c3339decSBarry Smith Logically Collective 540291d0ed5SPieter Ghysels 541291d0ed5SPieter Ghysels Input Parameters: 54229e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 543291d0ed5SPieter Ghysels 54429e0a805SPieter Ghysels Output Parameter: 54529e0a805SPieter Ghysels . leaf_size - Size of diagonal blocks in rank-structured approximation 546291d0ed5SPieter Ghysels 54729e0a805SPieter Ghysels Level: intermediate 548291d0ed5SPieter Ghysels 549291d0ed5SPieter Ghysels References: 55029e0a805SPieter Ghysels . * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/ 551291d0ed5SPieter Ghysels 55229e0a805SPieter Ghysels .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatSTRUMPACKSetCompLeafSize()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()` 553291d0ed5SPieter Ghysels @*/ 55429e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetCompLeafSize(Mat F, PetscInt *leaf_size) 555d71ae5a4SJacob Faibussowitsch { 556291d0ed5SPieter Ghysels PetscFunctionBegin; 557291d0ed5SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 55829e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKGetCompLeafSize_C", (Mat, PetscInt *), (F, leaf_size)); 55929e0a805SPieter Ghysels PetscValidLogicalCollectiveInt(F, *leaf_size, 2); 56029e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 56129e0a805SPieter Ghysels } 56229e0a805SPieter Ghysels 56329e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetGeometricNxyz_STRUMPACK(Mat F, PetscInt nx, PetscInt ny, PetscInt nz) 56429e0a805SPieter Ghysels { 56529e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 56629e0a805SPieter Ghysels 56729e0a805SPieter Ghysels PetscFunctionBegin; 56829e0a805SPieter Ghysels if (nx < 1) { 56929e0a805SPieter Ghysels PetscCheck(nx == PETSC_DECIDE || nx == PETSC_DEFAULT, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_OUTOFRANGE, "nx < 1"); 57029e0a805SPieter Ghysels nx = 1; 57129e0a805SPieter Ghysels } 57229e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_nx", STRUMPACK_set_nx(*S, nx)); 57329e0a805SPieter Ghysels if (ny < 1) { 57429e0a805SPieter Ghysels PetscCheck(ny == PETSC_DECIDE || ny == PETSC_DEFAULT, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_OUTOFRANGE, "ny < 1"); 57529e0a805SPieter Ghysels ny = 1; 57629e0a805SPieter Ghysels } 57729e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_ny", STRUMPACK_set_ny(*S, ny)); 57829e0a805SPieter Ghysels if (nz < 1) { 57929e0a805SPieter Ghysels PetscCheck(nz == PETSC_DECIDE || nz == PETSC_DEFAULT, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_OUTOFRANGE, "nz < 1"); 58029e0a805SPieter Ghysels nz = 1; 58129e0a805SPieter Ghysels } 58229e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_nz", STRUMPACK_set_nz(*S, nz)); 58329e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 58429e0a805SPieter Ghysels } 58529e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetGeometricComponents_STRUMPACK(Mat F, PetscInt nc) 58629e0a805SPieter Ghysels { 58729e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 58829e0a805SPieter Ghysels 58929e0a805SPieter Ghysels PetscFunctionBegin; 59029e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_components", STRUMPACK_set_components(*S, nc)); 59129e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 59229e0a805SPieter Ghysels } 59329e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetGeometricWidth_STRUMPACK(Mat F, PetscInt w) 59429e0a805SPieter Ghysels { 59529e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 59629e0a805SPieter Ghysels 59729e0a805SPieter Ghysels PetscFunctionBegin; 59829e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_separator_width", STRUMPACK_set_separator_width(*S, w)); 59929e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 60029e0a805SPieter Ghysels } 60129e0a805SPieter Ghysels 60229e0a805SPieter Ghysels /*@ 60329e0a805SPieter Ghysels MatSTRUMPACKSetGeometricNxyz - Set STRUMPACK mesh x, y and z dimensions, for use with GEOMETRIC ordering. 60429e0a805SPieter Ghysels 60529e0a805SPieter Ghysels Logically Collective 60629e0a805SPieter Ghysels 60729e0a805SPieter Ghysels Input Parameters: 60829e0a805SPieter Ghysels + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 60929e0a805SPieter Ghysels . nx - x dimension of the mesh 61029e0a805SPieter Ghysels . ny - y dimension of the mesh 61129e0a805SPieter Ghysels - nz - z dimension of the mesh 61229e0a805SPieter Ghysels 61329e0a805SPieter Ghysels Level: intermediate 61429e0a805SPieter Ghysels 61529e0a805SPieter Ghysels Notes: 61629e0a805SPieter Ghysels If the mesh is two (or one) dimensional one can use 1, `PETSC_DECIDE` or `PETSC_DEFAULT` 61729e0a805SPieter Ghysels for the missing z (and y) dimensions. 61829e0a805SPieter Ghysels 61929e0a805SPieter Ghysels References: 62029e0a805SPieter Ghysels . * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/ 62129e0a805SPieter Ghysels 62229e0a805SPieter Ghysels .seealso: `MatGetFactor()` 62329e0a805SPieter Ghysels @*/ 62429e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetGeometricNxyz(Mat F, PetscInt nx, PetscInt ny, PetscInt nz) 62529e0a805SPieter Ghysels { 62629e0a805SPieter Ghysels PetscFunctionBegin; 62729e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 62829e0a805SPieter Ghysels PetscValidLogicalCollectiveInt(F, nx, 2); 62929e0a805SPieter Ghysels PetscValidLogicalCollectiveInt(F, ny, 3); 63029e0a805SPieter Ghysels PetscValidLogicalCollectiveInt(F, nz, 4); 63129e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKSetGeometricNxyz_C", (Mat, PetscInt, PetscInt, PetscInt), (F, nx, ny, nz)); 63229e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 63329e0a805SPieter Ghysels } 63429e0a805SPieter Ghysels /*@ 63529e0a805SPieter Ghysels MatSTRUMPACKSetGeometricComponents - Set STRUMPACK number of degrees of freedom per mesh point, for use with GEOMETRIC ordering. 63629e0a805SPieter Ghysels 63729e0a805SPieter Ghysels Logically Collective 63829e0a805SPieter Ghysels 63929e0a805SPieter Ghysels Input Parameters: 64029e0a805SPieter Ghysels + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 64129e0a805SPieter Ghysels - nc - Number of components/dof's per grid point 64229e0a805SPieter Ghysels 64329e0a805SPieter Ghysels Options Database Key: 64429e0a805SPieter Ghysels . -mat_strumpack_geometric_components <1> - Number of components per mesh point, for geometric nested dissection ordering 64529e0a805SPieter Ghysels 64629e0a805SPieter Ghysels Level: intermediate 64729e0a805SPieter Ghysels 64829e0a805SPieter Ghysels References: 64929e0a805SPieter Ghysels . * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/ 65029e0a805SPieter Ghysels 65129e0a805SPieter Ghysels .seealso: `MatGetFactor()` 65229e0a805SPieter Ghysels @*/ 65329e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetGeometricComponents(Mat F, PetscInt nc) 65429e0a805SPieter Ghysels { 65529e0a805SPieter Ghysels PetscFunctionBegin; 65629e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 65729e0a805SPieter Ghysels PetscValidLogicalCollectiveInt(F, nc, 2); 65829e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKSetGeometricComponents_C", (Mat, PetscInt), (F, nc)); 65929e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 66029e0a805SPieter Ghysels } 66129e0a805SPieter Ghysels /*@ 66229e0a805SPieter Ghysels MatSTRUMPACKSetGeometricWidth - Set STRUMPACK width of the separator, for use with GEOMETRIC ordering. 66329e0a805SPieter Ghysels 66429e0a805SPieter Ghysels Logically Collective 66529e0a805SPieter Ghysels 66629e0a805SPieter Ghysels Input Parameters: 66729e0a805SPieter Ghysels + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 66829e0a805SPieter Ghysels - w - width of the separator 66929e0a805SPieter Ghysels 67029e0a805SPieter Ghysels Options Database Key: 67129e0a805SPieter Ghysels . -mat_strumpack_geometric_width <1> - Width of the separator of the mesh, for geometric nested dissection ordering 67229e0a805SPieter Ghysels 67329e0a805SPieter Ghysels Level: intermediate 67429e0a805SPieter Ghysels 67529e0a805SPieter Ghysels References: 67629e0a805SPieter Ghysels . * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/ 67729e0a805SPieter Ghysels 67829e0a805SPieter Ghysels .seealso: `MatGetFactor()` 67929e0a805SPieter Ghysels @*/ 68029e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetGeometricWidth(Mat F, PetscInt w) 68129e0a805SPieter Ghysels { 68229e0a805SPieter Ghysels PetscFunctionBegin; 68329e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 68429e0a805SPieter Ghysels PetscValidLogicalCollectiveInt(F, w, 2); 68529e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKSetGeometricWidth_C", (Mat, PetscInt), (F, w)); 68629e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 68729e0a805SPieter Ghysels } 68829e0a805SPieter Ghysels 68929e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetCompMinSepSize_STRUMPACK(Mat F, PetscInt min_sep_size) 69029e0a805SPieter Ghysels { 69129e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 69229e0a805SPieter Ghysels 69329e0a805SPieter Ghysels PetscFunctionBegin; 69429e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_compression_min_sep_size", STRUMPACK_set_compression_min_sep_size(*S, min_sep_size)); 69529e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 69629e0a805SPieter Ghysels } 69729e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetCompMinSepSize_STRUMPACK(Mat F, PetscInt *min_sep_size) 69829e0a805SPieter Ghysels { 69929e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 70029e0a805SPieter Ghysels 70129e0a805SPieter Ghysels PetscFunctionBegin; 70229e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression_min_sep_size", *min_sep_size = (PetscInt)STRUMPACK_compression_min_sep_size(*S)); 70329e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 70429e0a805SPieter Ghysels } 70529e0a805SPieter Ghysels 70629e0a805SPieter Ghysels /*@ 70729e0a805SPieter Ghysels MatSTRUMPACKSetCompMinSepSize - Set STRUMPACK minimum separator size for low-rank approximation 70829e0a805SPieter Ghysels 70929e0a805SPieter Ghysels Logically Collective 71029e0a805SPieter Ghysels 71129e0a805SPieter Ghysels Input Parameters: 71229e0a805SPieter Ghysels + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 71329e0a805SPieter Ghysels - min_sep_size - minimum dense matrix size for low-rank approximation 71429e0a805SPieter Ghysels 71529e0a805SPieter Ghysels Options Database Key: 71629e0a805SPieter Ghysels . -mat_strumpack_compression_min_sep_size <min_sep_size> - Minimum size of dense sub-block for low-rank compression 71729e0a805SPieter Ghysels 71829e0a805SPieter Ghysels Level: intermediate 71929e0a805SPieter Ghysels 72029e0a805SPieter Ghysels References: 72129e0a805SPieter Ghysels . * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/ 72229e0a805SPieter Ghysels 72329e0a805SPieter Ghysels .seealso: `MatGetFactor()`, `MatSTRUMPACKGetCompMinSepSize()` 72429e0a805SPieter Ghysels @*/ 72529e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetCompMinSepSize(Mat F, PetscInt min_sep_size) 72629e0a805SPieter Ghysels { 72729e0a805SPieter Ghysels PetscFunctionBegin; 72829e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 72929e0a805SPieter Ghysels PetscValidLogicalCollectiveInt(F, min_sep_size, 2); 73029e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKSetCompMinSepSize_C", (Mat, PetscInt), (F, min_sep_size)); 73129e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 73229e0a805SPieter Ghysels } 73329e0a805SPieter Ghysels /*@ 73429e0a805SPieter Ghysels MatSTRUMPACKGetCompMinSepSize - Get STRUMPACK minimum separator size for low-rank approximation 73529e0a805SPieter Ghysels 73629e0a805SPieter Ghysels Logically Collective 73729e0a805SPieter Ghysels 73829e0a805SPieter Ghysels Input Parameters: 73929e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 74029e0a805SPieter Ghysels 74129e0a805SPieter Ghysels Output Parameter: 74229e0a805SPieter Ghysels . min_sep_size - minimum dense matrix size for low-rank approximation 74329e0a805SPieter Ghysels 74429e0a805SPieter Ghysels Level: intermediate 74529e0a805SPieter Ghysels 74629e0a805SPieter Ghysels References: 74729e0a805SPieter Ghysels . * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/ 74829e0a805SPieter Ghysels 74929e0a805SPieter Ghysels .seealso: `MatGetFactor()`, `MatSTRUMPACKSetCompMinSepSize()` 75029e0a805SPieter Ghysels @*/ 75129e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetCompMinSepSize(Mat F, PetscInt *min_sep_size) 75229e0a805SPieter Ghysels { 75329e0a805SPieter Ghysels PetscFunctionBegin; 75429e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 75529e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKGetCompMinSepSize_C", (Mat, PetscInt *), (F, min_sep_size)); 75629e0a805SPieter Ghysels PetscValidLogicalCollectiveInt(F, *min_sep_size, 2); 75729e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 75829e0a805SPieter Ghysels } 75929e0a805SPieter Ghysels 76029e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetCompLossyPrecision_STRUMPACK(Mat F, PetscInt lossy_prec) 76129e0a805SPieter Ghysels { 76229e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 76329e0a805SPieter Ghysels 76429e0a805SPieter Ghysels PetscFunctionBegin; 76529e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_compression_lossy_precision", STRUMPACK_set_compression_lossy_precision(*S, lossy_prec)); 76629e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 76729e0a805SPieter Ghysels } 76829e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetCompLossyPrecision_STRUMPACK(Mat F, PetscInt *lossy_prec) 76929e0a805SPieter Ghysels { 77029e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 77129e0a805SPieter Ghysels 77229e0a805SPieter Ghysels PetscFunctionBegin; 77329e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression_lossy_precision", *lossy_prec = (PetscInt)STRUMPACK_compression_lossy_precision(*S)); 77429e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 77529e0a805SPieter Ghysels } 77629e0a805SPieter Ghysels 77729e0a805SPieter Ghysels /*@ 77829e0a805SPieter Ghysels MatSTRUMPACKSetCompLossyPrecision - Set STRUMPACK precision for lossy compression (requires ZFP support) 77929e0a805SPieter Ghysels 78029e0a805SPieter Ghysels Logically Collective 78129e0a805SPieter Ghysels 78229e0a805SPieter Ghysels Input Parameters: 78329e0a805SPieter Ghysels + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 78429e0a805SPieter Ghysels - lossy_prec - Number of bitplanes to use in lossy compression 78529e0a805SPieter Ghysels 78629e0a805SPieter Ghysels Options Database Key: 78729e0a805SPieter 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` 78829e0a805SPieter Ghysels 78929e0a805SPieter Ghysels Level: intermediate 79029e0a805SPieter Ghysels 79129e0a805SPieter Ghysels References: 79229e0a805SPieter Ghysels . * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/ 79329e0a805SPieter Ghysels 79429e0a805SPieter Ghysels .seealso: `MatGetFactor()`, `MatSTRUMPACKGetCompLossyPrecision()` 79529e0a805SPieter Ghysels @*/ 79629e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetCompLossyPrecision(Mat F, PetscInt lossy_prec) 79729e0a805SPieter Ghysels { 79829e0a805SPieter Ghysels PetscFunctionBegin; 79929e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 80029e0a805SPieter Ghysels PetscValidLogicalCollectiveInt(F, lossy_prec, 2); 80129e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKSetCompLossyPrecision_C", (Mat, PetscInt), (F, lossy_prec)); 80229e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 80329e0a805SPieter Ghysels } 80429e0a805SPieter Ghysels /*@ 80529e0a805SPieter Ghysels MatSTRUMPACKGetCompLossyPrecision - Get STRUMPACK precision for lossy compression (requires ZFP support) 80629e0a805SPieter Ghysels 80729e0a805SPieter Ghysels Logically Collective 80829e0a805SPieter Ghysels 80929e0a805SPieter Ghysels Input Parameters: 81029e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 81129e0a805SPieter Ghysels 81229e0a805SPieter Ghysels Output Parameter: 81329e0a805SPieter Ghysels . lossy_prec - Number of bitplanes to use in lossy compression 81429e0a805SPieter Ghysels 81529e0a805SPieter Ghysels Level: intermediate 81629e0a805SPieter Ghysels 81729e0a805SPieter Ghysels References: 81829e0a805SPieter Ghysels . * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/ 81929e0a805SPieter Ghysels 82029e0a805SPieter Ghysels .seealso: `MatGetFactor()`, `MatSTRUMPACKSetCompLossyPrecision()` 82129e0a805SPieter Ghysels @*/ 82229e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetCompLossyPrecision(Mat F, PetscInt *lossy_prec) 82329e0a805SPieter Ghysels { 82429e0a805SPieter Ghysels PetscFunctionBegin; 82529e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 82629e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKGetCompLossyPrecision_C", (Mat, PetscInt *), (F, lossy_prec)); 82729e0a805SPieter Ghysels PetscValidLogicalCollectiveInt(F, *lossy_prec, 2); 82829e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 82929e0a805SPieter Ghysels } 83029e0a805SPieter Ghysels 83129e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetCompButterflyLevels_STRUMPACK(Mat F, PetscInt bfly_lvls) 83229e0a805SPieter Ghysels { 83329e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 83429e0a805SPieter Ghysels 83529e0a805SPieter Ghysels PetscFunctionBegin; 83629e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_compression_butterfly_levels", STRUMPACK_set_compression_butterfly_levels(*S, bfly_lvls)); 83729e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 83829e0a805SPieter Ghysels } 83929e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetCompButterflyLevels_STRUMPACK(Mat F, PetscInt *bfly_lvls) 84029e0a805SPieter Ghysels { 84129e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 84229e0a805SPieter Ghysels 84329e0a805SPieter Ghysels PetscFunctionBegin; 84429e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression_butterfly_levels", *bfly_lvls = (PetscInt)STRUMPACK_compression_butterfly_levels(*S)); 84529e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 84629e0a805SPieter Ghysels } 84729e0a805SPieter Ghysels 84829e0a805SPieter Ghysels /*@ 84929e0a805SPieter Ghysels MatSTRUMPACKSetCompButterflyLevels - Set STRUMPACK number of butterfly levels in HODLR compression (requires ButterflyPACK support) 85029e0a805SPieter Ghysels 85129e0a805SPieter Ghysels Logically Collective 85229e0a805SPieter Ghysels 85329e0a805SPieter Ghysels Input Parameters: 85429e0a805SPieter Ghysels + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 85529e0a805SPieter Ghysels - bfly_lvls - Number of levels of butterfly compression in HODLR compression 85629e0a805SPieter Ghysels 85729e0a805SPieter Ghysels Options Database Key: 85829e0a805SPieter Ghysels . -mat_strumpack_compression_butterfly_levels <bfly_lvls> - Number of levels in the hierarchically off-diagonal matrix for which to use butterfly, when using `-pctype ilu`, (BLR_)HODLR compression 85929e0a805SPieter Ghysels 86029e0a805SPieter Ghysels Level: intermediate 86129e0a805SPieter Ghysels 86229e0a805SPieter Ghysels References: 86329e0a805SPieter Ghysels . * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/ 86429e0a805SPieter Ghysels 86529e0a805SPieter Ghysels .seealso: `MatGetFactor()`, `MatSTRUMPACKGetCompButterflyLevels()` 86629e0a805SPieter Ghysels @*/ 86729e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetCompButterflyLevels(Mat F, PetscInt bfly_lvls) 86829e0a805SPieter Ghysels { 86929e0a805SPieter Ghysels PetscFunctionBegin; 87029e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 87129e0a805SPieter Ghysels PetscValidLogicalCollectiveInt(F, bfly_lvls, 2); 87229e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKSetButterflyLevels_C", (Mat, PetscInt), (F, bfly_lvls)); 87329e0a805SPieter Ghysels PetscFunctionReturn(PETSC_SUCCESS); 87429e0a805SPieter Ghysels } 87529e0a805SPieter Ghysels /*@ 87629e0a805SPieter Ghysels MatSTRUMPACKGetCompButterflyLevels - Get STRUMPACK number of butterfly levels in HODLR compression (requires ButterflyPACK support) 87729e0a805SPieter Ghysels 87829e0a805SPieter Ghysels Logically Collective 87929e0a805SPieter Ghysels 88029e0a805SPieter Ghysels Input Parameters: 88129e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface 88229e0a805SPieter Ghysels 88329e0a805SPieter Ghysels Output Parameter: 88429e0a805SPieter Ghysels . bfly_lvls - Number of levels of butterfly compression in HODLR compression 88529e0a805SPieter Ghysels 88629e0a805SPieter Ghysels Level: intermediate 88729e0a805SPieter Ghysels 88829e0a805SPieter Ghysels References: 88929e0a805SPieter Ghysels . * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/ 89029e0a805SPieter Ghysels 89129e0a805SPieter Ghysels .seealso: `MatGetFactor()`, `MatSTRUMPACKSetCompButterflyLevels()` 89229e0a805SPieter Ghysels @*/ 89329e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetCompButterflyLevels(Mat F, PetscInt *bfly_lvls) 89429e0a805SPieter Ghysels { 89529e0a805SPieter Ghysels PetscFunctionBegin; 89629e0a805SPieter Ghysels PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 89729e0a805SPieter Ghysels PetscTryMethod(F, "MatSTRUMPACKGetButterflyLevels_C", (Mat, PetscInt *), (F, bfly_lvls)); 89829e0a805SPieter Ghysels PetscValidLogicalCollectiveInt(F, *bfly_lvls, 2); 8993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 900291d0ed5SPieter Ghysels } 901291d0ed5SPieter Ghysels 902d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolve_STRUMPACK(Mat A, Vec b_mpi, Vec x) 903d71ae5a4SJacob Faibussowitsch { 90435e8bcc9SJunchao Zhang STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)A->data; 9057d6ea485SPieter Ghysels STRUMPACK_RETURN_CODE sp_err; 9067d6ea485SPieter Ghysels const PetscScalar *bptr; 9077d6ea485SPieter Ghysels PetscScalar *xptr; 9087d6ea485SPieter Ghysels 9097d6ea485SPieter Ghysels PetscFunctionBegin; 9109566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &xptr)); 9119566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(b_mpi, &bptr)); 9120d08a34dSPieter Ghysels 913e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_solve", sp_err = STRUMPACK_solve(*S, (PetscScalar *)bptr, xptr, 0)); 9140d08a34dSPieter Ghysels switch (sp_err) { 915d71ae5a4SJacob Faibussowitsch case STRUMPACK_SUCCESS: 916d71ae5a4SJacob Faibussowitsch break; 9179371c9d4SSatish Balay case STRUMPACK_MATRIX_NOT_SET: { 9189371c9d4SSatish Balay SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix was not set"); 9199371c9d4SSatish Balay break; 9209371c9d4SSatish Balay } 9219371c9d4SSatish Balay case STRUMPACK_REORDERING_ERROR: { 9229371c9d4SSatish Balay SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix reordering failed"); 9239371c9d4SSatish Balay break; 9249371c9d4SSatish Balay } 925d71ae5a4SJacob Faibussowitsch default: 926d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: solve failed"); 9277d6ea485SPieter Ghysels } 9289566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &xptr)); 9299566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(b_mpi, &bptr)); 9303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9317d6ea485SPieter Ghysels } 9327d6ea485SPieter Ghysels 933d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMatSolve_STRUMPACK(Mat A, Mat B_mpi, Mat X) 934d71ae5a4SJacob Faibussowitsch { 93529e0a805SPieter Ghysels STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)A->data; 93629e0a805SPieter Ghysels STRUMPACK_RETURN_CODE sp_err; 9377d6ea485SPieter Ghysels PetscBool flg; 93829e0a805SPieter Ghysels PetscInt m = A->rmap->n, nrhs; 93929e0a805SPieter Ghysels const PetscScalar *bptr; 94029e0a805SPieter Ghysels PetscScalar *xptr; 9417d6ea485SPieter Ghysels 9427d6ea485SPieter Ghysels PetscFunctionBegin; 9439566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)B_mpi, &flg, MATSEQDENSE, MATMPIDENSE, NULL)); 9445f80ce2aSJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Matrix B must be MATDENSE matrix"); 9459566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)X, &flg, MATSEQDENSE, MATMPIDENSE, NULL)); 9465f80ce2aSJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Matrix X must be MATDENSE matrix"); 94729e0a805SPieter Ghysels 94829e0a805SPieter Ghysels PetscCall(MatGetSize(B_mpi, NULL, &nrhs)); 94929e0a805SPieter Ghysels PetscCall(MatDenseGetArray(X, &xptr)); 95029e0a805SPieter Ghysels PetscCall(MatDenseGetArrayRead(B_mpi, &bptr)); 95129e0a805SPieter Ghysels 95229e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_matsolve", sp_err = STRUMPACK_matsolve(*S, nrhs, bptr, m, xptr, m, 0)); 95329e0a805SPieter Ghysels switch (sp_err) { 95429e0a805SPieter Ghysels case STRUMPACK_SUCCESS: 95529e0a805SPieter Ghysels break; 95629e0a805SPieter Ghysels case STRUMPACK_MATRIX_NOT_SET: { 95729e0a805SPieter Ghysels SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix was not set"); 95829e0a805SPieter Ghysels break; 95929e0a805SPieter Ghysels } 96029e0a805SPieter Ghysels case STRUMPACK_REORDERING_ERROR: { 96129e0a805SPieter Ghysels SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix reordering failed"); 96229e0a805SPieter Ghysels break; 96329e0a805SPieter Ghysels } 96429e0a805SPieter Ghysels default: 96529e0a805SPieter Ghysels SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: solve failed"); 96629e0a805SPieter Ghysels } 96729e0a805SPieter Ghysels PetscCall(MatDenseRestoreArrayRead(B_mpi, &bptr)); 96829e0a805SPieter Ghysels PetscCall(MatDenseRestoreArray(X, &xptr)); 96929e0a805SPieter Ghysels 9703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9717d6ea485SPieter Ghysels } 9727d6ea485SPieter Ghysels 973d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_Info_STRUMPACK(Mat A, PetscViewer viewer) 974d71ae5a4SJacob Faibussowitsch { 975ad0c5e61SPieter Ghysels PetscFunctionBegin; 976ad0c5e61SPieter Ghysels /* check if matrix is strumpack type */ 9773ba16761SJacob Faibussowitsch if (A->ops->solve != MatSolve_STRUMPACK) PetscFunctionReturn(PETSC_SUCCESS); 9789566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "STRUMPACK sparse solver!\n")); 9793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 980ad0c5e61SPieter Ghysels } 981ad0c5e61SPieter Ghysels 982d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_STRUMPACK(Mat A, PetscViewer viewer) 983d71ae5a4SJacob Faibussowitsch { 984ad0c5e61SPieter Ghysels PetscBool iascii; 985ad0c5e61SPieter Ghysels PetscViewerFormat format; 986ad0c5e61SPieter Ghysels 987ad0c5e61SPieter Ghysels PetscFunctionBegin; 9889566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 989ad0c5e61SPieter Ghysels if (iascii) { 9909566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 9919566063dSJacob Faibussowitsch if (format == PETSC_VIEWER_ASCII_INFO) PetscCall(MatView_Info_STRUMPACK(A, viewer)); 992ad0c5e61SPieter Ghysels } 9933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 994ad0c5e61SPieter Ghysels } 9957d6ea485SPieter Ghysels 996d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLUFactorNumeric_STRUMPACK(Mat F, Mat A, const MatFactorInfo *info) 997d71ae5a4SJacob Faibussowitsch { 99835e8bcc9SJunchao Zhang STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data; 9997d6ea485SPieter Ghysels STRUMPACK_RETURN_CODE sp_err; 100029e0a805SPieter Ghysels Mat Aloc; 100129e0a805SPieter Ghysels const PetscScalar *av; 100229e0a805SPieter Ghysels const PetscInt *ai = NULL, *aj = NULL; 100329e0a805SPieter Ghysels PetscInt M = A->rmap->N, m = A->rmap->n, dummy; 100429e0a805SPieter Ghysels PetscBool ismpiaij, isseqaij, flg; 10057d6ea485SPieter Ghysels 10067d6ea485SPieter Ghysels PetscFunctionBegin; 100729e0a805SPieter Ghysels PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATSEQAIJ, &isseqaij)); 100829e0a805SPieter Ghysels PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATMPIAIJ, &ismpiaij)); 100929e0a805SPieter Ghysels if (ismpiaij) { 101029e0a805SPieter Ghysels PetscCall(MatMPIAIJGetLocalMat(A, MAT_INITIAL_MATRIX, &Aloc)); 101129e0a805SPieter Ghysels } else if (isseqaij) { 101229e0a805SPieter Ghysels PetscCall(PetscObjectReference((PetscObject)A)); 101329e0a805SPieter Ghysels Aloc = A; 101429e0a805SPieter Ghysels } else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not for type %s", ((PetscObject)A)->type_name); 101529e0a805SPieter Ghysels 101629e0a805SPieter Ghysels PetscCall(MatGetRowIJ(Aloc, 0, PETSC_FALSE, PETSC_FALSE, &dummy, &ai, &aj, &flg)); 101729e0a805SPieter Ghysels PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "GetRowIJ failed"); 101829e0a805SPieter Ghysels PetscCall(MatSeqAIJGetArrayRead(Aloc, &av)); 101929e0a805SPieter Ghysels 102029e0a805SPieter Ghysels if (ismpiaij) { 102129e0a805SPieter Ghysels const PetscInt *dist = NULL; 102229e0a805SPieter Ghysels PetscCall(MatGetOwnershipRanges(A, &dist)); 102329e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_distributed_csr_matrix", STRUMPACK_set_distributed_csr_matrix(*S, &m, ai, aj, av, dist, 0)); 102429e0a805SPieter Ghysels } else if (isseqaij) { 102529e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_csr_matrix", STRUMPACK_set_csr_matrix(*S, &M, ai, aj, av, 0)); 10267d6ea485SPieter Ghysels } 10277d6ea485SPieter Ghysels 102829e0a805SPieter Ghysels PetscCall(MatRestoreRowIJ(Aloc, 0, PETSC_FALSE, PETSC_FALSE, &dummy, &ai, &aj, &flg)); 102929e0a805SPieter Ghysels PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "RestoreRowIJ failed"); 103029e0a805SPieter Ghysels PetscCall(MatSeqAIJRestoreArrayRead(Aloc, &av)); 103129e0a805SPieter Ghysels PetscCall(MatDestroy(&Aloc)); 103229e0a805SPieter Ghysels 10337d6ea485SPieter Ghysels /* Reorder and Factor the matrix. */ 10347d6ea485SPieter Ghysels /* TODO figure out how to avoid reorder if the matrix values changed, but the pattern remains the same. */ 1035e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_reorder", sp_err = STRUMPACK_reorder(*S)); 1036e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_factor", sp_err = STRUMPACK_factor(*S)); 10370d08a34dSPieter Ghysels switch (sp_err) { 1038d71ae5a4SJacob Faibussowitsch case STRUMPACK_SUCCESS: 1039d71ae5a4SJacob Faibussowitsch break; 10409371c9d4SSatish Balay case STRUMPACK_MATRIX_NOT_SET: { 10419371c9d4SSatish Balay SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix was not set"); 10429371c9d4SSatish Balay break; 10439371c9d4SSatish Balay } 10449371c9d4SSatish Balay case STRUMPACK_REORDERING_ERROR: { 10459371c9d4SSatish Balay SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix reordering failed"); 10469371c9d4SSatish Balay break; 10479371c9d4SSatish Balay } 1048d71ae5a4SJacob Faibussowitsch default: 1049d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: factorization failed"); 10507d6ea485SPieter Ghysels } 1051cb250fa3SPieter Ghysels F->assembled = PETSC_TRUE; 1052cb250fa3SPieter Ghysels F->preallocated = PETSC_TRUE; 10533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10547d6ea485SPieter Ghysels } 10557d6ea485SPieter Ghysels 1056d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLUFactorSymbolic_STRUMPACK(Mat F, Mat A, IS r, IS c, const MatFactorInfo *info) 1057d71ae5a4SJacob Faibussowitsch { 10587d6ea485SPieter Ghysels PetscFunctionBegin; 10597d6ea485SPieter Ghysels F->ops->lufactornumeric = MatLUFactorNumeric_STRUMPACK; 10607d6ea485SPieter Ghysels F->ops->solve = MatSolve_STRUMPACK; 10617d6ea485SPieter Ghysels F->ops->matsolve = MatMatSolve_STRUMPACK; 10623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10637d6ea485SPieter Ghysels } 10647d6ea485SPieter Ghysels 1065d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorGetSolverType_aij_strumpack(Mat A, MatSolverType *type) 1066d71ae5a4SJacob Faibussowitsch { 10677d6ea485SPieter Ghysels PetscFunctionBegin; 10687d6ea485SPieter Ghysels *type = MATSOLVERSTRUMPACK; 10693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10707d6ea485SPieter Ghysels } 10717d6ea485SPieter Ghysels 1072575704cbSPieter Ghysels /*MC 107329e0a805SPieter Ghysels MATSOLVERSTRUMPACK = "strumpack" - A solver package providing a direct sparse solver (PCLU) 1074575704cbSPieter Ghysels and a preconditioner (PCILU) using low-rank compression via the external package STRUMPACK. 1075575704cbSPieter Ghysels 107629e0a805SPieter Ghysels Consult the STRUMPACK manual for more info, 107729e0a805SPieter Ghysels https://portal.nersc.gov/project/sparse/strumpack/master/ 1078575704cbSPieter Ghysels 107929e0a805SPieter Ghysels Use `./configure --download-strumpack --download-metis` to have PETSc installed with STRUMPACK. 108029e0a805SPieter Ghysels 108129e0a805SPieter Ghysels For full functionality, add `--download-slate --download-magma --download-parmetis --download-ptscotch --download-zfp --download-butterflypack`. 108229e0a805SPieter Ghysels SLATE provides GPU support in the multi-GPU setting, providing ScaLAPACK functionality but with GPU acceleration. 108329e0a805SPieter Ghysels MAGMA can optionally be used for on node GPU support instead cuBLAS/cuSOLVER, and performs slightly better. 108429e0a805SPieter Ghysels ParMETIS and PTScotch can be used for parallel fill-reducing ordering. 108529e0a805SPieter Ghysels ZFP is used for floating point compression of the sparse factors (LOSSY or LOSSLESS compression). 108629e0a805SPieter Ghysels ButterflyPACK is used for HODLR (Hierarchically Off-Diagonal Low Rank) and HODBF (Hierarchically Off-Diagonal Butterfly) compression of the sparse factors. 108729e0a805SPieter Ghysels 108829e0a805SPieter Ghysels Recommended use is 1 MPI rank per GPU. 1089575704cbSPieter Ghysels 10902ef1f0ffSBarry Smith Use `-pc_type lu` `-pc_factor_mat_solver_type strumpack` to use this as an exact (direct) solver. 1091575704cbSPieter Ghysels 109229e0a805SPieter Ghysels 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). 10932ef1f0ffSBarry Smith 10942ef1f0ffSBarry Smith Works with `MATAIJ` matrices 1095575704cbSPieter Ghysels 1096575704cbSPieter Ghysels Options Database Keys: 109729e0a805SPieter Ghysels + -mat_strumpack_verbose - Enable verbose output 109829e0a805SPieter 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 109929e0a805SPieter Ghysels . -mat_strumpack_compression_rel_tol - Relative compression tolerance, when using `-pctype ilu` 110029e0a805SPieter Ghysels . -mat_strumpack_compression_abs_tol - Absolute compression tolerance, when using `-pctype ilu` 110129e0a805SPieter Ghysels . -mat_strumpack_compression_min_sep_size - Minimum size of separator for rank-structured compression, when using `-pctype ilu` 110229e0a805SPieter Ghysels . -mat_strumpack_compression_leaf_size - Size of diagonal blocks in rank-structured approximation, when using `-pctype ilu` 110329e0a805SPieter Ghysels . -mat_strumpack_compression_lossy_precision - Precision when using lossy compression [1-64], when using `-pctype ilu`, compression LOSSY (requires ZFP support) 110429e0a805SPieter 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) 110529e0a805SPieter Ghysels . -mat_strumpack_gpu - Enable GPU acceleration in numerical factorization (not supported for all compression types) 110629e0a805SPieter Ghysels . -mat_strumpack_colperm <TRUE> - Permute matrix to make diagonal nonzeros 110729e0a805SPieter Ghysels . -mat_strumpack_reordering <METIS> - Sparsity reducing matrix reordering (choose one of) NATURAL METIS PARMETIS SCOTCH PTSCOTCH RCM GEOMETRIC AMD MMD AND MLF SPECTRAL 110829e0a805SPieter Ghysels . -mat_strumpack_geometric_xyz <1,1,1> - Mesh x,y,z dimensions, for use with GEOMETRIC ordering 110929e0a805SPieter Ghysels . -mat_strumpack_geometric_components <1> - Number of components per mesh point, for geometric nested dissection ordering 111029e0a805SPieter Ghysels . -mat_strumpack_geometric_width <1> - Width of the separator of the mesh, for geometric nested dissection ordering 111129e0a805SPieter Ghysels - -mat_strumpack_metis_nodeNDP - Use METIS_NodeNDP instead of METIS_NodeND, for a more balanced tree 1112575704cbSPieter Ghysels 1113575704cbSPieter Ghysels Level: beginner 1114575704cbSPieter Ghysels 111529e0a805SPieter Ghysels HODLR, BLR_HODBF and ZFP_BLR_HODLR compression require STRUMPACK to be configured with ButterflyPACK support (`--download-butterflypack`). 111629e0a805SPieter Ghysels 111729e0a805SPieter Ghysels LOSSY, LOSSLESS and ZFP_BLR_HODLR compression require STRUMPACK to be configured with ZFP support (`--download-zfp`). 111829e0a805SPieter Ghysels 11191cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `PCLU`, `PCILU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`, `PCFactorSetMatSolverType()`, `MatSolverType`, 112029e0a805SPieter Ghysels `MatGetFactor()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKReordering`, `MatSTRUMPACKCompressionType`, `MatSTRUMPACKSetColPerm()`. 1121575704cbSPieter Ghysels M*/ 1122d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetFactor_aij_strumpack(Mat A, MatFactorType ftype, Mat *F) 1123d71ae5a4SJacob Faibussowitsch { 11247d6ea485SPieter Ghysels Mat B; 11257d6ea485SPieter Ghysels PetscInt M = A->rmap->N, N = A->cmap->N; 112629e0a805SPieter Ghysels PetscBool verb, flg, set; 112729e0a805SPieter Ghysels PetscReal ctol; 112829e0a805SPieter Ghysels PetscInt min_sep_size, leaf_size, nxyz[3], nrdims, nc, w; 112929e0a805SPieter Ghysels #if defined(STRUMPACK_USE_ZFP) 113029e0a805SPieter Ghysels PetscInt lossy_prec; 113129e0a805SPieter Ghysels #endif 113229e0a805SPieter Ghysels #if defined(STRUMPACK_USE_BPACK) 113329e0a805SPieter Ghysels PetscInt bfly_lvls; 113429e0a805SPieter Ghysels #endif 113529e0a805SPieter Ghysels #if defined(STRUMPACK_USE_SLATE_SCALAPACK) 113629e0a805SPieter Ghysels PetscMPIInt mpithreads; 113729e0a805SPieter Ghysels #endif 113834f43fa5SPieter Ghysels STRUMPACK_SparseSolver *S; 113934f43fa5SPieter Ghysels STRUMPACK_INTERFACE iface; 114029e0a805SPieter Ghysels STRUMPACK_REORDERING_STRATEGY ndcurrent, ndvalue; 114129e0a805SPieter Ghysels STRUMPACK_COMPRESSION_TYPE compcurrent, compvalue; 11429371c9d4SSatish Balay const STRUMPACK_PRECISION table[2][2][2] = { 11439371c9d4SSatish Balay {{STRUMPACK_FLOATCOMPLEX_64, STRUMPACK_DOUBLECOMPLEX_64}, {STRUMPACK_FLOAT_64, STRUMPACK_DOUBLE_64}}, 11449371c9d4SSatish Balay {{STRUMPACK_FLOATCOMPLEX, STRUMPACK_DOUBLECOMPLEX}, {STRUMPACK_FLOAT, STRUMPACK_DOUBLE} } 11459371c9d4SSatish Balay }; 11464ac6704cSBarry Smith const STRUMPACK_PRECISION prec = table[(sizeof(PetscInt) == 8) ? 0 : 1][(PETSC_SCALAR == PETSC_COMPLEX) ? 0 : 1][(PETSC_REAL == PETSC_FLOAT) ? 0 : 1]; 114729e0a805SPieter Ghysels const char *const STRUMPACKNDTypes[] = {"NATURAL", "METIS", "PARMETIS", "SCOTCH", "PTSCOTCH", "RCM", "GEOMETRIC", "AMD", "MMD", "AND", "MLF", "SPECTRAL", "STRUMPACKNDTypes", "", 0}; 114829e0a805SPieter Ghysels const char *const CompTypes[] = {"NONE", "HSS", "BLR", "HODLR", "BLR_HODLR", "ZFP_BLR_HODLR", "LOSSLESS", "LOSSY", "CompTypes", "", 0}; 11497d6ea485SPieter Ghysels 11507d6ea485SPieter Ghysels PetscFunctionBegin; 115129e0a805SPieter Ghysels #if defined(STRUMPACK_USE_SLATE_SCALAPACK) 115229e0a805SPieter Ghysels PetscCallMPI(MPI_Query_thread(&mpithreads)); 115329e0a805SPieter Ghysels PetscCheck(mpithreads == MPI_THREAD_MULTIPLE, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP_SYS, "SLATE requires MPI_THREAD_MULTIPLE"); 115429e0a805SPieter Ghysels #endif 11557d6ea485SPieter Ghysels /* Create the factorization matrix */ 11569566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B)); 11579566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, M, N)); 115835e8bcc9SJunchao Zhang PetscCall(PetscStrallocpy("strumpack", &((PetscObject)B)->type_name)); 115935e8bcc9SJunchao Zhang PetscCall(MatSetUp(B)); 116029e0a805SPieter Ghysels PetscCall(MatSeqAIJSetPreallocation(B, 0, NULL)); 116129e0a805SPieter Ghysels PetscCall(MatMPIAIJSetPreallocation(B, 0, NULL, 0, NULL)); 116266e17bc3SBarry Smith B->trivialsymbolic = PETSC_TRUE; 1163575704cbSPieter Ghysels if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) { 11647d6ea485SPieter Ghysels B->ops->lufactorsymbolic = MatLUFactorSymbolic_STRUMPACK; 1165575704cbSPieter Ghysels B->ops->ilufactorsymbolic = MatLUFactorSymbolic_STRUMPACK; 1166575704cbSPieter Ghysels } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Factor type not supported"); 116735e8bcc9SJunchao Zhang B->ops->getinfo = MatGetInfo_External; 11687d6ea485SPieter Ghysels B->ops->view = MatView_STRUMPACK; 11697d6ea485SPieter Ghysels B->ops->destroy = MatDestroy_STRUMPACK; 11707d6ea485SPieter Ghysels B->ops->getdiagonal = MatGetDiagonal_STRUMPACK; 11719566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_aij_strumpack)); 11729566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetReordering_C", MatSTRUMPACKSetReordering_STRUMPACK)); 117329e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetReordering_C", MatSTRUMPACKGetReordering_STRUMPACK)); 11749566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetColPerm_C", MatSTRUMPACKSetColPerm_STRUMPACK)); 117529e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetColPerm_C", MatSTRUMPACKGetColPerm_STRUMPACK)); 117629e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetGeometricNxyz_C", MatSTRUMPACKSetGeometricNxyz_STRUMPACK)); 117729e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetGeometricComponents_C", MatSTRUMPACKSetGeometricComponents_STRUMPACK)); 117829e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetGeometricWidth_C", MatSTRUMPACKSetGeometricWidth_STRUMPACK)); 117929e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetGPU_C", MatSTRUMPACKSetGPU_STRUMPACK)); 118029e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetGPU_C", MatSTRUMPACKGetGPU_STRUMPACK)); 118129e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetCompression_C", MatSTRUMPACKSetCompression_STRUMPACK)); 118229e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetCompression_C", MatSTRUMPACKGetCompression_STRUMPACK)); 118329e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetCompRelTol_C", MatSTRUMPACKSetCompRelTol_STRUMPACK)); 118429e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetCompRelTol_C", MatSTRUMPACKGetCompRelTol_STRUMPACK)); 118529e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetCompAbsTol_C", MatSTRUMPACKSetCompAbsTol_STRUMPACK)); 118629e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetCompAbsTol_C", MatSTRUMPACKGetCompAbsTol_STRUMPACK)); 118729e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetCompLeafSize_C", MatSTRUMPACKSetCompLeafSize_STRUMPACK)); 118829e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetCompLeafSize_C", MatSTRUMPACKGetCompLeafSize_STRUMPACK)); 118929e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetCompMinSepSize_C", MatSTRUMPACKSetCompMinSepSize_STRUMPACK)); 119029e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetCompMinSepSize_C", MatSTRUMPACKGetCompMinSepSize_STRUMPACK)); 119129e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetCompLossyPrecision_C", MatSTRUMPACKSetCompLossyPrecision_STRUMPACK)); 119229e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetCompLossyPrecision_C", MatSTRUMPACKGetCompLossyPrecision_STRUMPACK)); 119329e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetCompButterflyLevels_C", MatSTRUMPACKSetCompButterflyLevels_STRUMPACK)); 119429e0a805SPieter Ghysels PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetCompButterflyLevels_C", MatSTRUMPACKGetCompButterflyLevels_STRUMPACK)); 1195575704cbSPieter Ghysels B->factortype = ftype; 119629e0a805SPieter Ghysels 119729e0a805SPieter Ghysels /* set solvertype */ 119829e0a805SPieter Ghysels PetscCall(PetscFree(B->solvertype)); 119929e0a805SPieter Ghysels PetscCall(PetscStrallocpy(MATSOLVERSTRUMPACK, &B->solvertype)); 120029e0a805SPieter Ghysels 12014dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&S)); 120235e8bcc9SJunchao Zhang B->data = S; 12030d08a34dSPieter Ghysels 120435e8bcc9SJunchao Zhang PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATSEQAIJ, &flg)); /* A might be MATSEQAIJCUSPARSE */ 12050d08a34dSPieter Ghysels iface = flg ? STRUMPACK_MT : STRUMPACK_MPI_DIST; 12067d6ea485SPieter Ghysels 120726cc229bSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)B), ((PetscObject)B)->prefix, "STRUMPACK Options", "Mat"); 1208575704cbSPieter Ghysels verb = PetscLogPrintInfo ? PETSC_TRUE : PETSC_FALSE; 12099566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_strumpack_verbose", "Print STRUMPACK information", "None", verb, &verb, NULL)); 12107d6ea485SPieter Ghysels 1211e77caa6dSBarry Smith PetscStackCallExternalVoid("STRUMPACK_init", STRUMPACK_init(S, PetscObjectComm((PetscObject)A), prec, iface, 0, NULL, verb)); 121255c022e5SPieter Ghysels 121329e0a805SPieter Ghysels /* By default, no compression is done. Compression is enabled when the user enables it with */ 121429e0a805SPieter Ghysels /* -mat_strumpack_compression with anything else than NONE, or when selecting ilu */ 121529e0a805SPieter Ghysels /* preconditioning, in which case we default to STRUMPACK_BLR compression. */ 121629e0a805SPieter Ghysels /* When compression is enabled, the STRUMPACK solver becomes an incomplete */ 121788382b05SPieter Ghysels /* (or approximate) LU factorization. */ 121829e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression", compcurrent = STRUMPACK_compression(*S)); 121929e0a805SPieter Ghysels PetscCall(PetscOptionsEnum("-mat_strumpack_compression", "Rank-structured compression type", "None", CompTypes, (PetscEnum)compcurrent, (PetscEnum *)&compvalue, &set)); 122029e0a805SPieter Ghysels if (set) { 122129e0a805SPieter Ghysels PetscCall(MatSTRUMPACKSetCompression(B, (MatSTRUMPACKCompressionType)compvalue)); 122229e0a805SPieter Ghysels } else { 122329e0a805SPieter Ghysels if (ftype == MAT_FACTOR_ILU) { PetscStackCallExternalVoid("STRUMPACK_set_compression", STRUMPACK_set_compression(*S, STRUMPACK_BLR)); } 122488382b05SPieter Ghysels } 122555c022e5SPieter Ghysels 122629e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression_rel_tol", ctol = (PetscReal)STRUMPACK_compression_rel_tol(*S)); 122729e0a805SPieter Ghysels PetscCall(PetscOptionsReal("-mat_strumpack_compression_rel_tol", "Relative compression tolerance", "None", ctol, &ctol, &set)); 122829e0a805SPieter Ghysels if (set) PetscStackCallExternalVoid("STRUMPACK_set_compression_rel_tol", STRUMPACK_set_compression_rel_tol(*S, (double)ctol)); 122929e0a805SPieter Ghysels 123029e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression_abs_tol", ctol = (PetscReal)STRUMPACK_compression_abs_tol(*S)); 123129e0a805SPieter Ghysels PetscCall(PetscOptionsReal("-mat_strumpack_compression_abs_tol", "Absolute compression tolerance", "None", ctol, &ctol, &set)); 123229e0a805SPieter Ghysels if (set) PetscStackCallExternalVoid("STRUMPACK_set_compression_abs_tol", STRUMPACK_set_compression_abs_tol(*S, (double)ctol)); 123329e0a805SPieter Ghysels 123429e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression_min_sep_size", min_sep_size = (PetscInt)STRUMPACK_compression_min_sep_size(*S)); 123529e0a805SPieter Ghysels PetscCall(PetscOptionsInt("-mat_strumpack_compression_min_sep_size", "Minimum size of separator for compression", "None", min_sep_size, &min_sep_size, &set)); 123629e0a805SPieter Ghysels if (set) PetscStackCallExternalVoid("STRUMPACK_set_compression_min_sep_size", STRUMPACK_set_compression_min_sep_size(*S, (int)min_sep_size)); 123729e0a805SPieter Ghysels 123829e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression_leaf_size", leaf_size = (PetscInt)STRUMPACK_compression_leaf_size(*S)); 123929e0a805SPieter Ghysels PetscCall(PetscOptionsInt("-mat_strumpack_compression_leaf_size", "Size of diagonal blocks in rank-structured approximation", "None", leaf_size, &leaf_size, &set)); 124029e0a805SPieter Ghysels if (set) PetscStackCallExternalVoid("STRUMPACK_set_compression_leaf_size", STRUMPACK_set_compression_leaf_size(*S, (int)leaf_size)); 124129e0a805SPieter Ghysels 124229e0a805SPieter Ghysels #if defined(STRUMPACK_USE_ZFP) 124329e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression_lossy_precision", lossy_prec = (PetscInt)STRUMPACK_compression_lossy_precision(*S)); 124429e0a805SPieter Ghysels PetscCall(PetscOptionsInt("-mat_strumpack_compression_lossy_precision", "Number of bitplanes to use in lossy compression", "None", lossy_prec, &lossy_prec, &set)); 124529e0a805SPieter Ghysels if (set) PetscStackCallExternalVoid("STRUMPACK_set_compression_lossy_precision", STRUMPACK_set_compression_lossy_precision(*S, (int)lossy_prec)); 124629e0a805SPieter Ghysels #endif 124729e0a805SPieter Ghysels 124829e0a805SPieter Ghysels #if defined(STRUMPACK_USE_BPACK) 124929e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_compression_butterfly_levels", bfly_lvls = (PetscInt)STRUMPACK_compression_butterfly_levels(*S)); 125029e0a805SPieter 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)); 125129e0a805SPieter Ghysels if (set) PetscStackCallExternalVoid("STRUMPACK_set_compression_butterfly_levels", STRUMPACK_set_compression_butterfly_levels(*S, (int)bfly_lvls)); 125229e0a805SPieter Ghysels #endif 125329e0a805SPieter Ghysels 125429e0a805SPieter Ghysels #if defined(STRUMPACK_USE_CUDA) || defined(STRUMPACK_USE_HIP) || defined(STRUMPACK_USE_SYCL) 125529e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_use_gpu", flg = (STRUMPACK_use_gpu(*S) == 0) ? PETSC_FALSE : PETSC_TRUE); 125629e0a805SPieter Ghysels PetscCall(PetscOptionsBool("-mat_strumpack_gpu", "Enable GPU acceleration (not supported for all compression types)", "None", flg, &flg, &set)); 125729e0a805SPieter Ghysels if (set) MatSTRUMPACKSetGPU(B, flg); 125829e0a805SPieter Ghysels #endif 125929e0a805SPieter Ghysels 126029e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_matching", flg = (STRUMPACK_matching(*S) == STRUMPACK_MATCHING_NONE) ? PETSC_FALSE : PETSC_TRUE); 126129e0a805SPieter Ghysels PetscCall(PetscOptionsBool("-mat_strumpack_colperm", "Find a col perm to get nonzero diagonal", "None", flg, &flg, &set)); 126229e0a805SPieter Ghysels if (set) PetscStackCallExternalVoid("STRUMPACK_set_matching", STRUMPACK_set_matching(*S, flg ? STRUMPACK_MATCHING_MAX_DIAGONAL_PRODUCT_SCALING : STRUMPACK_MATCHING_NONE)); 126329e0a805SPieter Ghysels 126429e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_reordering_method", ndcurrent = STRUMPACK_reordering_method(*S)); 126529e0a805SPieter Ghysels PetscCall(PetscOptionsEnum("-mat_strumpack_reordering", "Sparsity reducing matrix reordering", "None", STRUMPACKNDTypes, (PetscEnum)ndcurrent, (PetscEnum *)&ndvalue, &set)); 126629e0a805SPieter Ghysels if (set) PetscStackCallExternalVoid("STRUMPACK_set_reordering_method", STRUMPACK_set_reordering_method(*S, ndvalue)); 126729e0a805SPieter Ghysels 126829e0a805SPieter Ghysels /* geometric ordering, for a regular 1D/2D/3D mesh in the natural ordering, */ 126929e0a805SPieter Ghysels /* with nc DOF's per gridpoint, and possibly a wider stencil */ 127029e0a805SPieter Ghysels nrdims = 3; 127129e0a805SPieter Ghysels nxyz[0] = nxyz[1] = nxyz[2] = 1; 127229e0a805SPieter Ghysels PetscCall(PetscOptionsIntArray("-mat_strumpack_geometric_xyz", "Mesh sizes nx,ny,nz (Use 1 for default)", "", nxyz, &nrdims, &set)); 127329e0a805SPieter Ghysels if (set) { 127429e0a805SPieter 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."); 127529e0a805SPieter Ghysels PetscCall(MatSTRUMPACKSetGeometricNxyz(B, (int)nxyz[0], (int)nxyz[1], (int)nxyz[2])); 127629e0a805SPieter Ghysels } 127729e0a805SPieter Ghysels PetscCall(PetscOptionsInt("-mat_strumpack_geometric_components", "Number of components per mesh point, for geometric nested dissection ordering", "None", 1, &nc, &set)); 127829e0a805SPieter Ghysels if (set) PetscStackCallExternalVoid("STRUMPACK_set_components", STRUMPACK_set_components(*S, (int)nc)); 127929e0a805SPieter 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)); 128029e0a805SPieter Ghysels if (set) PetscStackCallExternalVoid("STRUMPACK_set_separator_width", STRUMPACK_set_separator_width(*S, (int)w)); 128129e0a805SPieter Ghysels 128229e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_use_METIS_NodeNDP", flg = (STRUMPACK_use_METIS_NodeNDP(*S) == 0) ? PETSC_FALSE : PETSC_TRUE); 128329e0a805SPieter Ghysels PetscCall(PetscOptionsBool("-mat_strumpack_metis_nodeNDP", "Use METIS_NodeNDP instead of METIS_NodeND, for a more balanced tree", "None", flg, &flg, &set)); 128429e0a805SPieter Ghysels if (set) { 128529e0a805SPieter Ghysels if (flg) { 128629e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_enable_METIS_NodeNDP", STRUMPACK_enable_METIS_NodeNDP(*S)); 128729e0a805SPieter Ghysels } else { 128829e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_disable_METIS_NodeNDP", STRUMPACK_disable_METIS_NodeNDP(*S)); 128929e0a805SPieter Ghysels } 129029e0a805SPieter Ghysels } 129129e0a805SPieter Ghysels 129229e0a805SPieter Ghysels /* Disable the outer iterative solver from STRUMPACK. */ 129329e0a805SPieter Ghysels /* When STRUMPACK is used as a direct solver, it will by default do iterative refinement. */ 129429e0a805SPieter Ghysels /* When STRUMPACK is used as an approximate factorization preconditioner (by enabling */ 129529e0a805SPieter Ghysels /* low-rank compression), it will use it's own preconditioned GMRES. Here we can disable */ 129629e0a805SPieter Ghysels /* the outer iterative solver, as PETSc uses STRUMPACK from within a KSP. */ 129729e0a805SPieter Ghysels PetscStackCallExternalVoid("STRUMPACK_set_Krylov_solver", STRUMPACK_set_Krylov_solver(*S, STRUMPACK_DIRECT)); 129829e0a805SPieter Ghysels 129929e0a805SPieter Ghysels PetscOptionsEnd(); 130026cc229bSBarry Smith 13017d6ea485SPieter Ghysels *F = B; 13023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13037d6ea485SPieter Ghysels } 13047d6ea485SPieter Ghysels 1305d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_STRUMPACK(void) 1306d71ae5a4SJacob Faibussowitsch { 13077d6ea485SPieter Ghysels PetscFunctionBegin; 13089566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSTRUMPACK, MATMPIAIJ, MAT_FACTOR_LU, MatGetFactor_aij_strumpack)); 13099566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSTRUMPACK, MATSEQAIJ, MAT_FACTOR_LU, MatGetFactor_aij_strumpack)); 13109566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSTRUMPACK, MATMPIAIJ, MAT_FACTOR_ILU, MatGetFactor_aij_strumpack)); 13119566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSTRUMPACK, MATSEQAIJ, MAT_FACTOR_ILU, MatGetFactor_aij_strumpack)); 13123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13137d6ea485SPieter Ghysels } 1314