xref: /petsc/src/mat/impls/aij/mpi/strumpack/strumpack.c (revision 29e0a8059c4f386829f8eba944f135cd960c2d7a)
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));
24*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKGetReordering_C", NULL));
259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetColPerm_C", NULL));
26*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKGetColPerm_C", NULL));
27*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetGeometricNxyz_C", NULL));
28*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetGeometricComponents_C", NULL));
29*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetGeometricWidth_C", NULL));
30*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetGPU_C", NULL));
31*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKGetGPU_C", NULL));
32*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetCompression_C", NULL));
33*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKGetCompression_C", NULL));
34*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetCompRelTol_C", NULL));
35*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKGetCompRelTol_C", NULL));
36*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetCompAbsTol_C", NULL));
37*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKGetCompAbsTol_C", NULL));
38*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetCompMaxRank_C", NULL));
39*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKGetCompMaxRank_C", NULL));
40*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetCompLeafSize_C", NULL));
41*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKGetCompLeafSize_C", NULL));
42*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetCompMinSepSize_C", NULL));
43*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKGetCompMinSepSize_C", NULL));
44*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetCompLossyPrecision_C", NULL));
45*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKGetCompLossyPrecision_C", NULL));
46*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetCompButterflyLevels_C", NULL));
47*29e0a805SPieter 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;
57*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_set_reordering_method", STRUMPACK_set_reordering_method(*S, (STRUMPACK_REORDERING_STRATEGY)reordering));
58*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
59*29e0a805SPieter Ghysels }
60*29e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetReordering_STRUMPACK(Mat F, MatSTRUMPACKReordering *reordering)
61*29e0a805SPieter Ghysels {
62*29e0a805SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
63*29e0a805SPieter Ghysels 
64*29e0a805SPieter Ghysels   PetscFunctionBegin;
65*29e0a805SPieter 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 
72*29e0a805SPieter Ghysels   Logically Collective
73*29e0a805SPieter Ghysels 
7434f43fa5SPieter Ghysels   Input Parameters:
75*29e0a805SPieter Ghysels + F          - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface
76*29e0a805SPieter 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 
81*29e0a805SPieter Ghysels   Level: intermediate
8234f43fa5SPieter Ghysels 
8334f43fa5SPieter Ghysels   References:
84*29e0a805SPieter Ghysels .  * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/
8534f43fa5SPieter Ghysels 
86*29e0a805SPieter 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 }
96*29e0a805SPieter Ghysels /*@
97*29e0a805SPieter Ghysels   MatSTRUMPACKGetReordering - Get STRUMPACK fill-reducing reordering
98*29e0a805SPieter Ghysels 
99*29e0a805SPieter Ghysels   Logically Collective
100*29e0a805SPieter Ghysels 
101*29e0a805SPieter Ghysels   Input Parameters:
102*29e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface
103*29e0a805SPieter Ghysels 
104*29e0a805SPieter Ghysels   Output Parameter:
105*29e0a805SPieter Ghysels . reordering - the code to be used to find the fill-reducing reordering
106*29e0a805SPieter Ghysels 
107*29e0a805SPieter Ghysels   Level: intermediate
108*29e0a805SPieter Ghysels 
109*29e0a805SPieter Ghysels   References:
110*29e0a805SPieter Ghysels .  * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/
111*29e0a805SPieter Ghysels 
112*29e0a805SPieter Ghysels .seealso: [](ch_matrices), `Mat`, `MatSTRUMPACKReordering`, `MatGetFactor()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()`
113*29e0a805SPieter Ghysels @*/
114*29e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetReordering(Mat F, MatSTRUMPACKReordering *reordering)
115*29e0a805SPieter Ghysels {
116*29e0a805SPieter Ghysels   PetscFunctionBegin;
117*29e0a805SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
118*29e0a805SPieter Ghysels   PetscTryMethod(F, "MatSTRUMPACKGetReordering_C", (Mat, MatSTRUMPACKReordering *), (F, reordering));
119*29e0a805SPieter Ghysels   PetscValidLogicalCollectiveEnum(F, *reordering, 2);
120*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
121*29e0a805SPieter 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;
128*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_set_matching", STRUMPACK_set_matching(*S, cperm ? STRUMPACK_MATCHING_MAX_DIAGONAL_PRODUCT_SCALING : STRUMPACK_MATCHING_NONE));
129*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
130*29e0a805SPieter Ghysels }
131*29e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetColPerm_STRUMPACK(Mat F, PetscBool *cperm)
132*29e0a805SPieter Ghysels {
133*29e0a805SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
134*29e0a805SPieter Ghysels 
135*29e0a805SPieter Ghysels   PetscFunctionBegin;
136*29e0a805SPieter 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 
152*29e0a805SPieter Ghysels   Level: intermediate
153575704cbSPieter Ghysels 
154575704cbSPieter Ghysels   References:
155*29e0a805SPieter Ghysels .  * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/
156575704cbSPieter Ghysels 
157*29e0a805SPieter 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 }
167*29e0a805SPieter Ghysels /*@
168*29e0a805SPieter Ghysels   MatSTRUMPACKGetColPerm - Get whether STRUMPACK will try to permute the columns of the matrix in order to get a nonzero diagonal
169575704cbSPieter Ghysels 
170*29e0a805SPieter Ghysels   Logically Collective
171*29e0a805SPieter Ghysels 
172*29e0a805SPieter Ghysels   Input Parameters:
173*29e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()`
174*29e0a805SPieter Ghysels 
175*29e0a805SPieter Ghysels   Output Parameter:
176*29e0a805SPieter Ghysels . cperm - Indicates whether STRUMPACK will permute columns
177*29e0a805SPieter Ghysels 
178*29e0a805SPieter Ghysels   Level: intermediate
179*29e0a805SPieter Ghysels 
180*29e0a805SPieter Ghysels   References:
181*29e0a805SPieter Ghysels .  * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/
182*29e0a805SPieter Ghysels 
183*29e0a805SPieter Ghysels .seealso: [](ch_matrices), `MatSTRUMPACKSetReordering()`, `Mat`, `MatGetFactor()`, `MatSTRUMPACKSetColPerm()`
184*29e0a805SPieter Ghysels @*/
185*29e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetColPerm(Mat F, PetscBool *cperm)
186*29e0a805SPieter Ghysels {
187*29e0a805SPieter Ghysels   PetscFunctionBegin;
188*29e0a805SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
189*29e0a805SPieter Ghysels   PetscTryMethod(F, "MatSTRUMPACKGetColPerm_C", (Mat, PetscBool *), (F, cperm));
190*29e0a805SPieter Ghysels   PetscValidLogicalCollectiveBool(F, *cperm, 2);
191*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
192*29e0a805SPieter Ghysels }
193*29e0a805SPieter Ghysels 
194*29e0a805SPieter 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;
199*29e0a805SPieter Ghysels   if (gpu) {
200*29e0a805SPieter Ghysels #if !(defined(STRUMPACK_USE_CUDA) || defined(STRUMPACK_USE_HIP) || defined(STRUMPACK_USE_SYCL))
201*29e0a805SPieter Ghysels     PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Warning: strumpack was not configured with GPU support\n"));
202*29e0a805SPieter Ghysels #endif
203*29e0a805SPieter Ghysels     PetscStackCallExternalVoid("STRUMPACK_enable_gpu", STRUMPACK_enable_gpu(*S));
204*29e0a805SPieter Ghysels   } else PetscStackCallExternalVoid("STRUMPACK_disable_gpu", STRUMPACK_disable_gpu(*S));
205*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
206*29e0a805SPieter Ghysels }
207*29e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetGPU_STRUMPACK(Mat F, PetscBool *gpu)
208*29e0a805SPieter Ghysels {
209*29e0a805SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
210*29e0a805SPieter Ghysels 
211*29e0a805SPieter Ghysels   PetscFunctionBegin;
212*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_use_gpu", *gpu = (PetscBool)STRUMPACK_use_gpu(*S));
2133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21434f43fa5SPieter Ghysels }
215e5a36eccSStefano Zampini 
21634f43fa5SPieter Ghysels /*@
217*29e0a805SPieter Ghysels   MatSTRUMPACKSetGPU - Set whether STRUMPACK should enable GPU acceleration (not supported for all compression types)
218*29e0a805SPieter Ghysels 
219*29e0a805SPieter Ghysels   Logically Collective
220*29e0a805SPieter Ghysels 
221*29e0a805SPieter Ghysels   Input Parameters:
222*29e0a805SPieter Ghysels + F   - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface
223*29e0a805SPieter Ghysels - gpu - whether or not to use GPU acceleration
224*29e0a805SPieter Ghysels 
225*29e0a805SPieter Ghysels   Options Database Key:
226*29e0a805SPieter Ghysels . -mat_strumpack_gpu <gpu> - true to use gpu offload
227*29e0a805SPieter Ghysels 
228*29e0a805SPieter Ghysels   Level: intermediate
229*29e0a805SPieter Ghysels 
230*29e0a805SPieter Ghysels   References:
231*29e0a805SPieter Ghysels .  * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/
232*29e0a805SPieter Ghysels 
233*29e0a805SPieter Ghysels .seealso: `MatGetFactor()`, `MatSTRUMPACKGetGPU()`
234*29e0a805SPieter Ghysels @*/
235*29e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetGPU(Mat F, PetscBool gpu)
236*29e0a805SPieter Ghysels {
237*29e0a805SPieter Ghysels   PetscFunctionBegin;
238*29e0a805SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
239*29e0a805SPieter Ghysels   PetscValidLogicalCollectiveBool(F, gpu, 2);
240*29e0a805SPieter Ghysels   PetscTryMethod(F, "MatSTRUMPACKSetGPU_C", (Mat, PetscBool), (F, gpu));
241*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
242*29e0a805SPieter Ghysels }
243*29e0a805SPieter Ghysels /*@
244*29e0a805SPieter Ghysels   MatSTRUMPACKGetGPU - Get whether STRUMPACK will try to use GPU acceleration (not supported for all compression types)
245*29e0a805SPieter Ghysels 
246*29e0a805SPieter Ghysels   Logically Collective
247*29e0a805SPieter Ghysels 
248*29e0a805SPieter Ghysels   Input Parameters:
249*29e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface
250*29e0a805SPieter Ghysels 
251*29e0a805SPieter Ghysels   Output Parameter:
252*29e0a805SPieter Ghysels . gpu - whether or not STRUMPACK will try to use GPU acceleration
253*29e0a805SPieter Ghysels 
254*29e0a805SPieter Ghysels   Level: intermediate
255*29e0a805SPieter Ghysels 
256*29e0a805SPieter Ghysels   References:
257*29e0a805SPieter Ghysels .  * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/
258*29e0a805SPieter Ghysels 
259*29e0a805SPieter Ghysels .seealso: `MatGetFactor()`, `MatSTRUMPACKSetGPU()`
260*29e0a805SPieter Ghysels @*/
261*29e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetGPU(Mat F, PetscBool *gpu)
262*29e0a805SPieter Ghysels {
263*29e0a805SPieter Ghysels   PetscFunctionBegin;
264*29e0a805SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
265*29e0a805SPieter Ghysels   PetscTryMethod(F, "MatSTRUMPACKGetGPU_C", (Mat, PetscBool *), (F, gpu));
266*29e0a805SPieter Ghysels   PetscValidLogicalCollectiveBool(F, *gpu, 2);
267*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
268*29e0a805SPieter Ghysels }
269*29e0a805SPieter Ghysels 
270*29e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetCompression_STRUMPACK(Mat F, MatSTRUMPACKCompressionType comp)
271*29e0a805SPieter Ghysels {
272*29e0a805SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
273*29e0a805SPieter Ghysels 
274*29e0a805SPieter Ghysels   PetscFunctionBegin;
275*29e0a805SPieter Ghysels #if !defined(STRUMPACK_HAVE_BPACK)
276*29e0a805SPieter 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");
277*29e0a805SPieter Ghysels #endif
278*29e0a805SPieter Ghysels #if !defined(STRUMPACK_HAVE_ZFP)
279*29e0a805SPieter 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");
280*29e0a805SPieter Ghysels #endif
281*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_set_compression", STRUMPACK_set_compression(*S, (STRUMPACK_COMPRESSION_TYPE)comp));
282*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
283*29e0a805SPieter Ghysels }
284*29e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetCompression_STRUMPACK(Mat F, MatSTRUMPACKCompressionType *comp)
285*29e0a805SPieter Ghysels {
286*29e0a805SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
287*29e0a805SPieter Ghysels 
288*29e0a805SPieter Ghysels   PetscFunctionBegin;
289*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_compression", *comp = (MatSTRUMPACKCompressionType)STRUMPACK_compression(*S));
290*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
291*29e0a805SPieter Ghysels }
292*29e0a805SPieter Ghysels 
293*29e0a805SPieter Ghysels /*@
294*29e0a805SPieter Ghysels   MatSTRUMPACKSetCompression - Set STRUMPACK compression type
295*29e0a805SPieter Ghysels 
296*29e0a805SPieter Ghysels   Input Parameters:
297*29e0a805SPieter Ghysels + F    - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface
298*29e0a805SPieter Ghysels - comp - Type of compression to be used in the approximate sparse factorization
299*29e0a805SPieter Ghysels 
300*29e0a805SPieter Ghysels   Options Database Key:
301*29e0a805SPieter 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
302*29e0a805SPieter Ghysels 
303*29e0a805SPieter Ghysels   Level: intermediate
304*29e0a805SPieter Ghysels 
305*29e0a805SPieter Ghysels   Notes:
306*29e0a805SPieter Ghysels   Default for `comp` is `MAT_STRUMPACK_COMPRESSION_TYPE_NONE` for `-pc_type lu` and `MAT_STRUMPACK_COMPRESSION_TYPE_BLR`
307*29e0a805SPieter Ghysels   for `-pc_type ilu`
308*29e0a805SPieter Ghysels 
309*29e0a805SPieter Ghysels   References:
310*29e0a805SPieter Ghysels .  * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/
311*29e0a805SPieter Ghysels 
312*29e0a805SPieter Ghysels .seealso: `MatGetFactor()`, `MatSTRUMPACKCompressionType`, `MatSTRUMPACKGetCompression()`
313*29e0a805SPieter Ghysels @*/
314*29e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetCompression(Mat F, MatSTRUMPACKCompressionType comp)
315*29e0a805SPieter Ghysels {
316*29e0a805SPieter Ghysels   PetscFunctionBegin;
317*29e0a805SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
318*29e0a805SPieter Ghysels   PetscValidLogicalCollectiveEnum(F, comp, 2);
319*29e0a805SPieter Ghysels   PetscTryMethod(F, "MatSTRUMPACKSetCompression_C", (Mat, MatSTRUMPACKCompressionType), (F, comp));
320*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
321*29e0a805SPieter Ghysels }
322*29e0a805SPieter Ghysels /*@
323*29e0a805SPieter Ghysels   MatSTRUMPACKGetCompression - Get STRUMPACK compression type
324*29e0a805SPieter Ghysels 
325*29e0a805SPieter Ghysels   Input Parameters:
326*29e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface
327*29e0a805SPieter Ghysels 
328*29e0a805SPieter Ghysels   Output Parameter:
329*29e0a805SPieter Ghysels . comp - Type of compression to be used in the approximate sparse factorization
330*29e0a805SPieter Ghysels 
331*29e0a805SPieter Ghysels   Level: intermediate
332*29e0a805SPieter Ghysels 
333*29e0a805SPieter Ghysels   Notes:
334*29e0a805SPieter Ghysels   Default is `MAT_STRUMPACK_COMPRESSION_TYPE_NONE` for `-pc_type lu` and `MAT_STRUMPACK_COMPRESSION_TYPE_BLR` for `-pc_type ilu`
335*29e0a805SPieter Ghysels 
336*29e0a805SPieter Ghysels   References:
337*29e0a805SPieter Ghysels .  * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/
338*29e0a805SPieter Ghysels 
339*29e0a805SPieter Ghysels .seealso: `MatGetFactor()`, `MatSTRUMPACKCompressionType`, `MatSTRUMPACKSetCompression()`
340*29e0a805SPieter Ghysels @*/
341*29e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetCompression(Mat F, MatSTRUMPACKCompressionType *comp)
342*29e0a805SPieter Ghysels {
343*29e0a805SPieter Ghysels   PetscFunctionBegin;
344*29e0a805SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
345*29e0a805SPieter Ghysels   PetscTryMethod(F, "MatSTRUMPACKGetCompression_C", (Mat, MatSTRUMPACKCompressionType *), (F, comp));
346*29e0a805SPieter Ghysels   PetscValidLogicalCollectiveEnum(F, *comp, 2);
347*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
348*29e0a805SPieter Ghysels }
349*29e0a805SPieter Ghysels 
350*29e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetCompRelTol_STRUMPACK(Mat F, PetscReal rtol)
351*29e0a805SPieter Ghysels {
352*29e0a805SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
353*29e0a805SPieter Ghysels 
354*29e0a805SPieter Ghysels   PetscFunctionBegin;
355*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_set_compression_rel_tol", STRUMPACK_set_compression_rel_tol(*S, rtol));
356*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
357*29e0a805SPieter Ghysels }
358*29e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetCompRelTol_STRUMPACK(Mat F, PetscReal *rtol)
359*29e0a805SPieter Ghysels {
360*29e0a805SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
361*29e0a805SPieter Ghysels 
362*29e0a805SPieter Ghysels   PetscFunctionBegin;
363*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_compression_rel_tol", *rtol = (PetscReal)STRUMPACK_compression_rel_tol(*S));
364*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
365*29e0a805SPieter Ghysels }
366*29e0a805SPieter Ghysels 
367*29e0a805SPieter Ghysels /*@
368*29e0a805SPieter 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:
377*29e0a805SPieter Ghysels . -mat_strumpack_compression_rel_tol <1e-4> - Relative compression tolerance, when using `-pctype ilu`
37834f43fa5SPieter Ghysels 
379*29e0a805SPieter Ghysels   Level: intermediate
38034f43fa5SPieter Ghysels 
38134f43fa5SPieter Ghysels   References:
382*29e0a805SPieter Ghysels .  * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/
38334f43fa5SPieter Ghysels 
384*29e0a805SPieter Ghysels .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatSTRUMPACKGetCompRelTol()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()`
38534f43fa5SPieter Ghysels @*/
386*29e0a805SPieter 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);
391*29e0a805SPieter Ghysels   PetscTryMethod(F, "MatSTRUMPACKSetCompRelTol_C", (Mat, PetscReal), (F, rtol));
392*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
393*29e0a805SPieter Ghysels }
394*29e0a805SPieter Ghysels /*@
395*29e0a805SPieter Ghysels   MatSTRUMPACKGetCompRelTol - Get STRUMPACK relative tolerance for compression
396*29e0a805SPieter Ghysels 
397*29e0a805SPieter Ghysels   Logically Collective
398*29e0a805SPieter Ghysels 
399*29e0a805SPieter Ghysels   Input Parameters:
400*29e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()`
401*29e0a805SPieter Ghysels 
402*29e0a805SPieter Ghysels   Output Parameter:
403*29e0a805SPieter Ghysels . rtol - relative compression tolerance
404*29e0a805SPieter Ghysels 
405*29e0a805SPieter Ghysels   Level: intermediate
406*29e0a805SPieter Ghysels 
407*29e0a805SPieter Ghysels   References:
408*29e0a805SPieter Ghysels .  * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/
409*29e0a805SPieter Ghysels 
410*29e0a805SPieter Ghysels .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatSTRUMPACKSetCompRelTol()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()`
411*29e0a805SPieter Ghysels @*/
412*29e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetCompRelTol(Mat F, PetscReal *rtol)
413*29e0a805SPieter Ghysels {
414*29e0a805SPieter Ghysels   PetscFunctionBegin;
415*29e0a805SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
416*29e0a805SPieter Ghysels   PetscTryMethod(F, "MatSTRUMPACKGetCompRelTol_C", (Mat, PetscReal *), (F, rtol));
417*29e0a805SPieter Ghysels   PetscValidLogicalCollectiveReal(F, *rtol, 2);
4183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41934f43fa5SPieter Ghysels }
42034f43fa5SPieter Ghysels 
421*29e0a805SPieter 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;
426*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_set_compression_abs_tol", STRUMPACK_set_compression_abs_tol(*S, atol));
427*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
428*29e0a805SPieter Ghysels }
429*29e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetCompAbsTol_STRUMPACK(Mat F, PetscReal *atol)
430*29e0a805SPieter Ghysels {
431*29e0a805SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
432*29e0a805SPieter Ghysels 
433*29e0a805SPieter Ghysels   PetscFunctionBegin;
434*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_compression_abs_tol", *atol = (PetscReal)STRUMPACK_compression_abs_tol(*S));
4353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
43634f43fa5SPieter Ghysels }
437e5a36eccSStefano Zampini 
43834f43fa5SPieter Ghysels /*@
439*29e0a805SPieter 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:
448*29e0a805SPieter Ghysels . -mat_strumpack_compression_abs_tol <1e-10> - Absolute compression tolerance, when using `-pctype ilu`
44934f43fa5SPieter Ghysels 
450*29e0a805SPieter Ghysels   Level: intermediate
45134f43fa5SPieter Ghysels 
45234f43fa5SPieter Ghysels   References:
453*29e0a805SPieter Ghysels .  * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/
45434f43fa5SPieter Ghysels 
455*29e0a805SPieter Ghysels .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatSTRUMPACKGetCompAbsTol()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()`
45634f43fa5SPieter Ghysels @*/
457*29e0a805SPieter 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);
462*29e0a805SPieter Ghysels   PetscTryMethod(F, "MatSTRUMPACKSetCompAbsTol_C", (Mat, PetscReal), (F, atol));
4633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
46434f43fa5SPieter Ghysels }
46534f43fa5SPieter Ghysels /*@
466*29e0a805SPieter Ghysels   MatSTRUMPACKGetCompAbsTol - Get STRUMPACK absolute tolerance for compression
467147403d9SBarry Smith 
468c3339decSBarry Smith   Logically Collective
46934f43fa5SPieter Ghysels 
47034f43fa5SPieter Ghysels   Input Parameters:
471*29e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()`
47234f43fa5SPieter Ghysels 
473*29e0a805SPieter Ghysels   Output Parameter:
474*29e0a805SPieter Ghysels . atol - absolute compression tolerance
47534f43fa5SPieter Ghysels 
476*29e0a805SPieter Ghysels   Level: intermediate
47734f43fa5SPieter Ghysels 
47834f43fa5SPieter Ghysels   References:
479*29e0a805SPieter Ghysels .  * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/
48034f43fa5SPieter Ghysels 
481*29e0a805SPieter Ghysels .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatSTRUMPACKSetCompAbsTol()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()`
48234f43fa5SPieter Ghysels @*/
483*29e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetCompAbsTol(Mat F, PetscReal *atol)
484d71ae5a4SJacob Faibussowitsch {
48534f43fa5SPieter Ghysels   PetscFunctionBegin;
48634f43fa5SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
487*29e0a805SPieter Ghysels   PetscTryMethod(F, "MatSTRUMPACKGetCompAbsTol_C", (Mat, PetscReal *), (F, atol));
488*29e0a805SPieter Ghysels   PetscValidLogicalCollectiveReal(F, *atol, 2);
4893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
49034f43fa5SPieter Ghysels }
49134f43fa5SPieter Ghysels 
492*29e0a805SPieter 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;
497*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_set_compression_leaf_size", STRUMPACK_set_compression_leaf_size(*S, leaf_size));
498*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
499*29e0a805SPieter Ghysels }
500*29e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetCompLeafSize_STRUMPACK(Mat F, PetscInt *leaf_size)
501*29e0a805SPieter Ghysels {
502*29e0a805SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
503*29e0a805SPieter Ghysels 
504*29e0a805SPieter Ghysels   PetscFunctionBegin;
505*29e0a805SPieter 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 /*@
510*29e0a805SPieter Ghysels   MatSTRUMPACKSetCompLeafSize - Set STRUMPACK leaf size for HSS, BLR, HODLR...
511147403d9SBarry Smith 
512c3339decSBarry Smith   Logically Collective
513a36bf211SPieter Ghysels 
514a36bf211SPieter Ghysels   Input Parameters:
515*29e0a805SPieter Ghysels + F         - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface
516*29e0a805SPieter Ghysels - leaf_size - Size of diagonal blocks in rank-structured approximation
517a36bf211SPieter Ghysels 
51811a5261eSBarry Smith   Options Database Key:
519*29e0a805SPieter Ghysels . -mat_strumpack_compression_leaf_size - Size of diagonal blocks in rank-structured approximation, when using `-pctype ilu`
520a36bf211SPieter Ghysels 
521*29e0a805SPieter Ghysels   Level: intermediate
522a36bf211SPieter Ghysels 
523a36bf211SPieter Ghysels   References:
524*29e0a805SPieter Ghysels .  * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/
525a36bf211SPieter Ghysels 
526*29e0a805SPieter Ghysels .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatSTRUMPACKGetCompLeafSize()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()`
527a36bf211SPieter Ghysels @*/
528*29e0a805SPieter 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);
533*29e0a805SPieter Ghysels   PetscTryMethod(F, "MatSTRUMPACKSetCompLeafSize_C", (Mat, PetscInt), (F, leaf_size));
5343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
535a36bf211SPieter Ghysels }
536291d0ed5SPieter Ghysels /*@
537*29e0a805SPieter Ghysels   MatSTRUMPACKGetCompLeafSize - Get STRUMPACK leaf size for HSS, BLR, HODLR...
538147403d9SBarry Smith 
539c3339decSBarry Smith   Logically Collective
540291d0ed5SPieter Ghysels 
541291d0ed5SPieter Ghysels   Input Parameters:
542*29e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface
543291d0ed5SPieter Ghysels 
544*29e0a805SPieter Ghysels   Output Parameter:
545*29e0a805SPieter Ghysels . leaf_size - Size of diagonal blocks in rank-structured approximation
546291d0ed5SPieter Ghysels 
547*29e0a805SPieter Ghysels   Level: intermediate
548291d0ed5SPieter Ghysels 
549291d0ed5SPieter Ghysels   References:
550*29e0a805SPieter Ghysels .  * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/
551291d0ed5SPieter Ghysels 
552*29e0a805SPieter Ghysels .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatSTRUMPACKSetCompLeafSize()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()`
553291d0ed5SPieter Ghysels @*/
554*29e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetCompLeafSize(Mat F, PetscInt *leaf_size)
555d71ae5a4SJacob Faibussowitsch {
556291d0ed5SPieter Ghysels   PetscFunctionBegin;
557291d0ed5SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
558*29e0a805SPieter Ghysels   PetscTryMethod(F, "MatSTRUMPACKGetCompLeafSize_C", (Mat, PetscInt *), (F, leaf_size));
559*29e0a805SPieter Ghysels   PetscValidLogicalCollectiveInt(F, *leaf_size, 2);
560*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
561*29e0a805SPieter Ghysels }
562*29e0a805SPieter Ghysels 
563*29e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetGeometricNxyz_STRUMPACK(Mat F, PetscInt nx, PetscInt ny, PetscInt nz)
564*29e0a805SPieter Ghysels {
565*29e0a805SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
566*29e0a805SPieter Ghysels 
567*29e0a805SPieter Ghysels   PetscFunctionBegin;
568*29e0a805SPieter Ghysels   if (nx < 1) {
569*29e0a805SPieter Ghysels     PetscCheck(nx == PETSC_DECIDE || nx == PETSC_DEFAULT, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_OUTOFRANGE, "nx < 1");
570*29e0a805SPieter Ghysels     nx = 1;
571*29e0a805SPieter Ghysels   }
572*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_set_nx", STRUMPACK_set_nx(*S, nx));
573*29e0a805SPieter Ghysels   if (ny < 1) {
574*29e0a805SPieter Ghysels     PetscCheck(ny == PETSC_DECIDE || ny == PETSC_DEFAULT, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_OUTOFRANGE, "ny < 1");
575*29e0a805SPieter Ghysels     ny = 1;
576*29e0a805SPieter Ghysels   }
577*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_set_ny", STRUMPACK_set_ny(*S, ny));
578*29e0a805SPieter Ghysels   if (nz < 1) {
579*29e0a805SPieter Ghysels     PetscCheck(nz == PETSC_DECIDE || nz == PETSC_DEFAULT, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_OUTOFRANGE, "nz < 1");
580*29e0a805SPieter Ghysels     nz = 1;
581*29e0a805SPieter Ghysels   }
582*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_set_nz", STRUMPACK_set_nz(*S, nz));
583*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
584*29e0a805SPieter Ghysels }
585*29e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetGeometricComponents_STRUMPACK(Mat F, PetscInt nc)
586*29e0a805SPieter Ghysels {
587*29e0a805SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
588*29e0a805SPieter Ghysels 
589*29e0a805SPieter Ghysels   PetscFunctionBegin;
590*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_set_components", STRUMPACK_set_components(*S, nc));
591*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
592*29e0a805SPieter Ghysels }
593*29e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetGeometricWidth_STRUMPACK(Mat F, PetscInt w)
594*29e0a805SPieter Ghysels {
595*29e0a805SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
596*29e0a805SPieter Ghysels 
597*29e0a805SPieter Ghysels   PetscFunctionBegin;
598*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_set_separator_width", STRUMPACK_set_separator_width(*S, w));
599*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
600*29e0a805SPieter Ghysels }
601*29e0a805SPieter Ghysels 
602*29e0a805SPieter Ghysels /*@
603*29e0a805SPieter Ghysels   MatSTRUMPACKSetGeometricNxyz - Set STRUMPACK mesh x, y and z dimensions, for use with GEOMETRIC ordering.
604*29e0a805SPieter Ghysels 
605*29e0a805SPieter Ghysels   Logically Collective
606*29e0a805SPieter Ghysels 
607*29e0a805SPieter Ghysels   Input Parameters:
608*29e0a805SPieter Ghysels + F  - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface
609*29e0a805SPieter Ghysels . nx - x dimension of the mesh
610*29e0a805SPieter Ghysels . ny - y dimension of the mesh
611*29e0a805SPieter Ghysels - nz - z dimension of the mesh
612*29e0a805SPieter Ghysels 
613*29e0a805SPieter Ghysels   Level: intermediate
614*29e0a805SPieter Ghysels 
615*29e0a805SPieter Ghysels   Notes:
616*29e0a805SPieter Ghysels   If the mesh is two (or one) dimensional one can use 1, `PETSC_DECIDE` or `PETSC_DEFAULT`
617*29e0a805SPieter Ghysels   for the missing z (and y) dimensions.
618*29e0a805SPieter Ghysels 
619*29e0a805SPieter Ghysels   References:
620*29e0a805SPieter Ghysels .  * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/
621*29e0a805SPieter Ghysels 
622*29e0a805SPieter Ghysels .seealso: `MatGetFactor()`
623*29e0a805SPieter Ghysels @*/
624*29e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetGeometricNxyz(Mat F, PetscInt nx, PetscInt ny, PetscInt nz)
625*29e0a805SPieter Ghysels {
626*29e0a805SPieter Ghysels   PetscFunctionBegin;
627*29e0a805SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
628*29e0a805SPieter Ghysels   PetscValidLogicalCollectiveInt(F, nx, 2);
629*29e0a805SPieter Ghysels   PetscValidLogicalCollectiveInt(F, ny, 3);
630*29e0a805SPieter Ghysels   PetscValidLogicalCollectiveInt(F, nz, 4);
631*29e0a805SPieter Ghysels   PetscTryMethod(F, "MatSTRUMPACKSetGeometricNxyz_C", (Mat, PetscInt, PetscInt, PetscInt), (F, nx, ny, nz));
632*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
633*29e0a805SPieter Ghysels }
634*29e0a805SPieter Ghysels /*@
635*29e0a805SPieter Ghysels   MatSTRUMPACKSetGeometricComponents - Set STRUMPACK number of degrees of freedom per mesh point, for use with GEOMETRIC ordering.
636*29e0a805SPieter Ghysels 
637*29e0a805SPieter Ghysels   Logically Collective
638*29e0a805SPieter Ghysels 
639*29e0a805SPieter Ghysels   Input Parameters:
640*29e0a805SPieter Ghysels + F  - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface
641*29e0a805SPieter Ghysels - nc - Number of components/dof's per grid point
642*29e0a805SPieter Ghysels 
643*29e0a805SPieter Ghysels   Options Database Key:
644*29e0a805SPieter Ghysels . -mat_strumpack_geometric_components <1> - Number of components per mesh point, for geometric nested dissection ordering
645*29e0a805SPieter Ghysels 
646*29e0a805SPieter Ghysels   Level: intermediate
647*29e0a805SPieter Ghysels 
648*29e0a805SPieter Ghysels   References:
649*29e0a805SPieter Ghysels .  * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/
650*29e0a805SPieter Ghysels 
651*29e0a805SPieter Ghysels .seealso: `MatGetFactor()`
652*29e0a805SPieter Ghysels @*/
653*29e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetGeometricComponents(Mat F, PetscInt nc)
654*29e0a805SPieter Ghysels {
655*29e0a805SPieter Ghysels   PetscFunctionBegin;
656*29e0a805SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
657*29e0a805SPieter Ghysels   PetscValidLogicalCollectiveInt(F, nc, 2);
658*29e0a805SPieter Ghysels   PetscTryMethod(F, "MatSTRUMPACKSetGeometricComponents_C", (Mat, PetscInt), (F, nc));
659*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
660*29e0a805SPieter Ghysels }
661*29e0a805SPieter Ghysels /*@
662*29e0a805SPieter Ghysels   MatSTRUMPACKSetGeometricWidth - Set STRUMPACK width of the separator, for use with GEOMETRIC ordering.
663*29e0a805SPieter Ghysels 
664*29e0a805SPieter Ghysels   Logically Collective
665*29e0a805SPieter Ghysels 
666*29e0a805SPieter Ghysels   Input Parameters:
667*29e0a805SPieter Ghysels + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface
668*29e0a805SPieter Ghysels - w - width of the separator
669*29e0a805SPieter Ghysels 
670*29e0a805SPieter Ghysels   Options Database Key:
671*29e0a805SPieter Ghysels . -mat_strumpack_geometric_width <1> - Width of the separator of the mesh, for geometric nested dissection ordering
672*29e0a805SPieter Ghysels 
673*29e0a805SPieter Ghysels   Level: intermediate
674*29e0a805SPieter Ghysels 
675*29e0a805SPieter Ghysels   References:
676*29e0a805SPieter Ghysels .  * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/
677*29e0a805SPieter Ghysels 
678*29e0a805SPieter Ghysels .seealso: `MatGetFactor()`
679*29e0a805SPieter Ghysels @*/
680*29e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetGeometricWidth(Mat F, PetscInt w)
681*29e0a805SPieter Ghysels {
682*29e0a805SPieter Ghysels   PetscFunctionBegin;
683*29e0a805SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
684*29e0a805SPieter Ghysels   PetscValidLogicalCollectiveInt(F, w, 2);
685*29e0a805SPieter Ghysels   PetscTryMethod(F, "MatSTRUMPACKSetGeometricWidth_C", (Mat, PetscInt), (F, w));
686*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
687*29e0a805SPieter Ghysels }
688*29e0a805SPieter Ghysels 
689*29e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetCompMinSepSize_STRUMPACK(Mat F, PetscInt min_sep_size)
690*29e0a805SPieter Ghysels {
691*29e0a805SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
692*29e0a805SPieter Ghysels 
693*29e0a805SPieter Ghysels   PetscFunctionBegin;
694*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_set_compression_min_sep_size", STRUMPACK_set_compression_min_sep_size(*S, min_sep_size));
695*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
696*29e0a805SPieter Ghysels }
697*29e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetCompMinSepSize_STRUMPACK(Mat F, PetscInt *min_sep_size)
698*29e0a805SPieter Ghysels {
699*29e0a805SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
700*29e0a805SPieter Ghysels 
701*29e0a805SPieter Ghysels   PetscFunctionBegin;
702*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_compression_min_sep_size", *min_sep_size = (PetscInt)STRUMPACK_compression_min_sep_size(*S));
703*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
704*29e0a805SPieter Ghysels }
705*29e0a805SPieter Ghysels 
706*29e0a805SPieter Ghysels /*@
707*29e0a805SPieter Ghysels   MatSTRUMPACKSetCompMinSepSize - Set STRUMPACK minimum separator size for low-rank approximation
708*29e0a805SPieter Ghysels 
709*29e0a805SPieter Ghysels   Logically Collective
710*29e0a805SPieter Ghysels 
711*29e0a805SPieter Ghysels   Input Parameters:
712*29e0a805SPieter Ghysels + F            - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface
713*29e0a805SPieter Ghysels - min_sep_size - minimum dense matrix size for low-rank approximation
714*29e0a805SPieter Ghysels 
715*29e0a805SPieter Ghysels   Options Database Key:
716*29e0a805SPieter Ghysels . -mat_strumpack_compression_min_sep_size <min_sep_size> - Minimum size of dense sub-block for low-rank compression
717*29e0a805SPieter Ghysels 
718*29e0a805SPieter Ghysels   Level: intermediate
719*29e0a805SPieter Ghysels 
720*29e0a805SPieter Ghysels   References:
721*29e0a805SPieter Ghysels .  * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/
722*29e0a805SPieter Ghysels 
723*29e0a805SPieter Ghysels .seealso: `MatGetFactor()`, `MatSTRUMPACKGetCompMinSepSize()`
724*29e0a805SPieter Ghysels @*/
725*29e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetCompMinSepSize(Mat F, PetscInt min_sep_size)
726*29e0a805SPieter Ghysels {
727*29e0a805SPieter Ghysels   PetscFunctionBegin;
728*29e0a805SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
729*29e0a805SPieter Ghysels   PetscValidLogicalCollectiveInt(F, min_sep_size, 2);
730*29e0a805SPieter Ghysels   PetscTryMethod(F, "MatSTRUMPACKSetCompMinSepSize_C", (Mat, PetscInt), (F, min_sep_size));
731*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
732*29e0a805SPieter Ghysels }
733*29e0a805SPieter Ghysels /*@
734*29e0a805SPieter Ghysels   MatSTRUMPACKGetCompMinSepSize - Get STRUMPACK minimum separator size for low-rank approximation
735*29e0a805SPieter Ghysels 
736*29e0a805SPieter Ghysels   Logically Collective
737*29e0a805SPieter Ghysels 
738*29e0a805SPieter Ghysels   Input Parameters:
739*29e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface
740*29e0a805SPieter Ghysels 
741*29e0a805SPieter Ghysels   Output Parameter:
742*29e0a805SPieter Ghysels . min_sep_size - minimum dense matrix size for low-rank approximation
743*29e0a805SPieter Ghysels 
744*29e0a805SPieter Ghysels   Level: intermediate
745*29e0a805SPieter Ghysels 
746*29e0a805SPieter Ghysels   References:
747*29e0a805SPieter Ghysels .  * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/
748*29e0a805SPieter Ghysels 
749*29e0a805SPieter Ghysels .seealso: `MatGetFactor()`, `MatSTRUMPACKSetCompMinSepSize()`
750*29e0a805SPieter Ghysels @*/
751*29e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetCompMinSepSize(Mat F, PetscInt *min_sep_size)
752*29e0a805SPieter Ghysels {
753*29e0a805SPieter Ghysels   PetscFunctionBegin;
754*29e0a805SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
755*29e0a805SPieter Ghysels   PetscTryMethod(F, "MatSTRUMPACKGetCompMinSepSize_C", (Mat, PetscInt *), (F, min_sep_size));
756*29e0a805SPieter Ghysels   PetscValidLogicalCollectiveInt(F, *min_sep_size, 2);
757*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
758*29e0a805SPieter Ghysels }
759*29e0a805SPieter Ghysels 
760*29e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetCompLossyPrecision_STRUMPACK(Mat F, PetscInt lossy_prec)
761*29e0a805SPieter Ghysels {
762*29e0a805SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
763*29e0a805SPieter Ghysels 
764*29e0a805SPieter Ghysels   PetscFunctionBegin;
765*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_set_compression_lossy_precision", STRUMPACK_set_compression_lossy_precision(*S, lossy_prec));
766*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
767*29e0a805SPieter Ghysels }
768*29e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetCompLossyPrecision_STRUMPACK(Mat F, PetscInt *lossy_prec)
769*29e0a805SPieter Ghysels {
770*29e0a805SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
771*29e0a805SPieter Ghysels 
772*29e0a805SPieter Ghysels   PetscFunctionBegin;
773*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_compression_lossy_precision", *lossy_prec = (PetscInt)STRUMPACK_compression_lossy_precision(*S));
774*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
775*29e0a805SPieter Ghysels }
776*29e0a805SPieter Ghysels 
777*29e0a805SPieter Ghysels /*@
778*29e0a805SPieter Ghysels   MatSTRUMPACKSetCompLossyPrecision - Set STRUMPACK precision for lossy compression (requires ZFP support)
779*29e0a805SPieter Ghysels 
780*29e0a805SPieter Ghysels   Logically Collective
781*29e0a805SPieter Ghysels 
782*29e0a805SPieter Ghysels   Input Parameters:
783*29e0a805SPieter Ghysels + F          - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface
784*29e0a805SPieter Ghysels - lossy_prec - Number of bitplanes to use in lossy compression
785*29e0a805SPieter Ghysels 
786*29e0a805SPieter Ghysels   Options Database Key:
787*29e0a805SPieter 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`
788*29e0a805SPieter Ghysels 
789*29e0a805SPieter Ghysels   Level: intermediate
790*29e0a805SPieter Ghysels 
791*29e0a805SPieter Ghysels   References:
792*29e0a805SPieter Ghysels .  * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/
793*29e0a805SPieter Ghysels 
794*29e0a805SPieter Ghysels .seealso: `MatGetFactor()`, `MatSTRUMPACKGetCompLossyPrecision()`
795*29e0a805SPieter Ghysels @*/
796*29e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetCompLossyPrecision(Mat F, PetscInt lossy_prec)
797*29e0a805SPieter Ghysels {
798*29e0a805SPieter Ghysels   PetscFunctionBegin;
799*29e0a805SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
800*29e0a805SPieter Ghysels   PetscValidLogicalCollectiveInt(F, lossy_prec, 2);
801*29e0a805SPieter Ghysels   PetscTryMethod(F, "MatSTRUMPACKSetCompLossyPrecision_C", (Mat, PetscInt), (F, lossy_prec));
802*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
803*29e0a805SPieter Ghysels }
804*29e0a805SPieter Ghysels /*@
805*29e0a805SPieter Ghysels   MatSTRUMPACKGetCompLossyPrecision - Get STRUMPACK precision for lossy compression (requires ZFP support)
806*29e0a805SPieter Ghysels 
807*29e0a805SPieter Ghysels   Logically Collective
808*29e0a805SPieter Ghysels 
809*29e0a805SPieter Ghysels   Input Parameters:
810*29e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface
811*29e0a805SPieter Ghysels 
812*29e0a805SPieter Ghysels   Output Parameter:
813*29e0a805SPieter Ghysels . lossy_prec - Number of bitplanes to use in lossy compression
814*29e0a805SPieter Ghysels 
815*29e0a805SPieter Ghysels   Level: intermediate
816*29e0a805SPieter Ghysels 
817*29e0a805SPieter Ghysels   References:
818*29e0a805SPieter Ghysels .  * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/
819*29e0a805SPieter Ghysels 
820*29e0a805SPieter Ghysels .seealso: `MatGetFactor()`, `MatSTRUMPACKSetCompLossyPrecision()`
821*29e0a805SPieter Ghysels @*/
822*29e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetCompLossyPrecision(Mat F, PetscInt *lossy_prec)
823*29e0a805SPieter Ghysels {
824*29e0a805SPieter Ghysels   PetscFunctionBegin;
825*29e0a805SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
826*29e0a805SPieter Ghysels   PetscTryMethod(F, "MatSTRUMPACKGetCompLossyPrecision_C", (Mat, PetscInt *), (F, lossy_prec));
827*29e0a805SPieter Ghysels   PetscValidLogicalCollectiveInt(F, *lossy_prec, 2);
828*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
829*29e0a805SPieter Ghysels }
830*29e0a805SPieter Ghysels 
831*29e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetCompButterflyLevels_STRUMPACK(Mat F, PetscInt bfly_lvls)
832*29e0a805SPieter Ghysels {
833*29e0a805SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
834*29e0a805SPieter Ghysels 
835*29e0a805SPieter Ghysels   PetscFunctionBegin;
836*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_set_compression_butterfly_levels", STRUMPACK_set_compression_butterfly_levels(*S, bfly_lvls));
837*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
838*29e0a805SPieter Ghysels }
839*29e0a805SPieter Ghysels static PetscErrorCode MatSTRUMPACKGetCompButterflyLevels_STRUMPACK(Mat F, PetscInt *bfly_lvls)
840*29e0a805SPieter Ghysels {
841*29e0a805SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
842*29e0a805SPieter Ghysels 
843*29e0a805SPieter Ghysels   PetscFunctionBegin;
844*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_compression_butterfly_levels", *bfly_lvls = (PetscInt)STRUMPACK_compression_butterfly_levels(*S));
845*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
846*29e0a805SPieter Ghysels }
847*29e0a805SPieter Ghysels 
848*29e0a805SPieter Ghysels /*@
849*29e0a805SPieter Ghysels   MatSTRUMPACKSetCompButterflyLevels - Set STRUMPACK number of butterfly levels in HODLR compression (requires ButterflyPACK support)
850*29e0a805SPieter Ghysels 
851*29e0a805SPieter Ghysels   Logically Collective
852*29e0a805SPieter Ghysels 
853*29e0a805SPieter Ghysels   Input Parameters:
854*29e0a805SPieter Ghysels + F         - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface
855*29e0a805SPieter Ghysels - bfly_lvls - Number of levels of butterfly compression in HODLR compression
856*29e0a805SPieter Ghysels 
857*29e0a805SPieter Ghysels   Options Database Key:
858*29e0a805SPieter 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
859*29e0a805SPieter Ghysels 
860*29e0a805SPieter Ghysels   Level: intermediate
861*29e0a805SPieter Ghysels 
862*29e0a805SPieter Ghysels   References:
863*29e0a805SPieter Ghysels .  * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/
864*29e0a805SPieter Ghysels 
865*29e0a805SPieter Ghysels .seealso: `MatGetFactor()`, `MatSTRUMPACKGetCompButterflyLevels()`
866*29e0a805SPieter Ghysels @*/
867*29e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKSetCompButterflyLevels(Mat F, PetscInt bfly_lvls)
868*29e0a805SPieter Ghysels {
869*29e0a805SPieter Ghysels   PetscFunctionBegin;
870*29e0a805SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
871*29e0a805SPieter Ghysels   PetscValidLogicalCollectiveInt(F, bfly_lvls, 2);
872*29e0a805SPieter Ghysels   PetscTryMethod(F, "MatSTRUMPACKSetButterflyLevels_C", (Mat, PetscInt), (F, bfly_lvls));
873*29e0a805SPieter Ghysels   PetscFunctionReturn(PETSC_SUCCESS);
874*29e0a805SPieter Ghysels }
875*29e0a805SPieter Ghysels /*@
876*29e0a805SPieter Ghysels   MatSTRUMPACKGetCompButterflyLevels - Get STRUMPACK number of butterfly levels in HODLR compression (requires ButterflyPACK support)
877*29e0a805SPieter Ghysels 
878*29e0a805SPieter Ghysels   Logically Collective
879*29e0a805SPieter Ghysels 
880*29e0a805SPieter Ghysels   Input Parameters:
881*29e0a805SPieter Ghysels . F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface
882*29e0a805SPieter Ghysels 
883*29e0a805SPieter Ghysels   Output Parameter:
884*29e0a805SPieter Ghysels . bfly_lvls - Number of levels of butterfly compression in HODLR compression
885*29e0a805SPieter Ghysels 
886*29e0a805SPieter Ghysels   Level: intermediate
887*29e0a805SPieter Ghysels 
888*29e0a805SPieter Ghysels   References:
889*29e0a805SPieter Ghysels .  * - STRUMPACK documentation: https://portal.nersc.gov/project/sparse/strumpack/master/
890*29e0a805SPieter Ghysels 
891*29e0a805SPieter Ghysels .seealso: `MatGetFactor()`, `MatSTRUMPACKSetCompButterflyLevels()`
892*29e0a805SPieter Ghysels @*/
893*29e0a805SPieter Ghysels PetscErrorCode MatSTRUMPACKGetCompButterflyLevels(Mat F, PetscInt *bfly_lvls)
894*29e0a805SPieter Ghysels {
895*29e0a805SPieter Ghysels   PetscFunctionBegin;
896*29e0a805SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
897*29e0a805SPieter Ghysels   PetscTryMethod(F, "MatSTRUMPACKGetButterflyLevels_C", (Mat, PetscInt *), (F, bfly_lvls));
898*29e0a805SPieter 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 {
935*29e0a805SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)A->data;
936*29e0a805SPieter Ghysels   STRUMPACK_RETURN_CODE   sp_err;
9377d6ea485SPieter Ghysels   PetscBool               flg;
938*29e0a805SPieter Ghysels   PetscInt                m = A->rmap->n, nrhs;
939*29e0a805SPieter Ghysels   const PetscScalar      *bptr;
940*29e0a805SPieter 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");
947*29e0a805SPieter Ghysels 
948*29e0a805SPieter Ghysels   PetscCall(MatGetSize(B_mpi, NULL, &nrhs));
949*29e0a805SPieter Ghysels   PetscCall(MatDenseGetArray(X, &xptr));
950*29e0a805SPieter Ghysels   PetscCall(MatDenseGetArrayRead(B_mpi, &bptr));
951*29e0a805SPieter Ghysels 
952*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_matsolve", sp_err = STRUMPACK_matsolve(*S, nrhs, bptr, m, xptr, m, 0));
953*29e0a805SPieter Ghysels   switch (sp_err) {
954*29e0a805SPieter Ghysels   case STRUMPACK_SUCCESS:
955*29e0a805SPieter Ghysels     break;
956*29e0a805SPieter Ghysels   case STRUMPACK_MATRIX_NOT_SET: {
957*29e0a805SPieter Ghysels     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix was not set");
958*29e0a805SPieter Ghysels     break;
959*29e0a805SPieter Ghysels   }
960*29e0a805SPieter Ghysels   case STRUMPACK_REORDERING_ERROR: {
961*29e0a805SPieter Ghysels     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix reordering failed");
962*29e0a805SPieter Ghysels     break;
963*29e0a805SPieter Ghysels   }
964*29e0a805SPieter Ghysels   default:
965*29e0a805SPieter Ghysels     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: solve failed");
966*29e0a805SPieter Ghysels   }
967*29e0a805SPieter Ghysels   PetscCall(MatDenseRestoreArrayRead(B_mpi, &bptr));
968*29e0a805SPieter Ghysels   PetscCall(MatDenseRestoreArray(X, &xptr));
969*29e0a805SPieter 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;
1000*29e0a805SPieter Ghysels   Mat                     Aloc;
1001*29e0a805SPieter Ghysels   const PetscScalar      *av;
1002*29e0a805SPieter Ghysels   const PetscInt         *ai = NULL, *aj = NULL;
1003*29e0a805SPieter Ghysels   PetscInt                M = A->rmap->N, m = A->rmap->n, dummy;
1004*29e0a805SPieter Ghysels   PetscBool               ismpiaij, isseqaij, flg;
10057d6ea485SPieter Ghysels 
10067d6ea485SPieter Ghysels   PetscFunctionBegin;
1007*29e0a805SPieter Ghysels   PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATSEQAIJ, &isseqaij));
1008*29e0a805SPieter Ghysels   PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATMPIAIJ, &ismpiaij));
1009*29e0a805SPieter Ghysels   if (ismpiaij) {
1010*29e0a805SPieter Ghysels     PetscCall(MatMPIAIJGetLocalMat(A, MAT_INITIAL_MATRIX, &Aloc));
1011*29e0a805SPieter Ghysels   } else if (isseqaij) {
1012*29e0a805SPieter Ghysels     PetscCall(PetscObjectReference((PetscObject)A));
1013*29e0a805SPieter Ghysels     Aloc = A;
1014*29e0a805SPieter Ghysels   } else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not for type %s", ((PetscObject)A)->type_name);
1015*29e0a805SPieter Ghysels 
1016*29e0a805SPieter Ghysels   PetscCall(MatGetRowIJ(Aloc, 0, PETSC_FALSE, PETSC_FALSE, &dummy, &ai, &aj, &flg));
1017*29e0a805SPieter Ghysels   PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "GetRowIJ failed");
1018*29e0a805SPieter Ghysels   PetscCall(MatSeqAIJGetArrayRead(Aloc, &av));
1019*29e0a805SPieter Ghysels 
1020*29e0a805SPieter Ghysels   if (ismpiaij) {
1021*29e0a805SPieter Ghysels     const PetscInt *dist = NULL;
1022*29e0a805SPieter Ghysels     PetscCall(MatGetOwnershipRanges(A, &dist));
1023*29e0a805SPieter Ghysels     PetscStackCallExternalVoid("STRUMPACK_set_distributed_csr_matrix", STRUMPACK_set_distributed_csr_matrix(*S, &m, ai, aj, av, dist, 0));
1024*29e0a805SPieter Ghysels   } else if (isseqaij) {
1025*29e0a805SPieter Ghysels     PetscStackCallExternalVoid("STRUMPACK_set_csr_matrix", STRUMPACK_set_csr_matrix(*S, &M, ai, aj, av, 0));
10267d6ea485SPieter Ghysels   }
10277d6ea485SPieter Ghysels 
1028*29e0a805SPieter Ghysels   PetscCall(MatRestoreRowIJ(Aloc, 0, PETSC_FALSE, PETSC_FALSE, &dummy, &ai, &aj, &flg));
1029*29e0a805SPieter Ghysels   PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "RestoreRowIJ failed");
1030*29e0a805SPieter Ghysels   PetscCall(MatSeqAIJRestoreArrayRead(Aloc, &av));
1031*29e0a805SPieter Ghysels   PetscCall(MatDestroy(&Aloc));
1032*29e0a805SPieter 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
1073*29e0a805SPieter 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 
1076*29e0a805SPieter Ghysels   Consult the STRUMPACK manual for more info,
1077*29e0a805SPieter Ghysels     https://portal.nersc.gov/project/sparse/strumpack/master/
1078575704cbSPieter Ghysels 
1079*29e0a805SPieter Ghysels   Use `./configure --download-strumpack --download-metis` to have PETSc installed with STRUMPACK.
1080*29e0a805SPieter Ghysels 
1081*29e0a805SPieter Ghysels   For full functionality, add `--download-slate --download-magma --download-parmetis --download-ptscotch --download-zfp --download-butterflypack`.
1082*29e0a805SPieter Ghysels   SLATE provides GPU support in the multi-GPU setting, providing ScaLAPACK functionality but with GPU acceleration.
1083*29e0a805SPieter Ghysels   MAGMA can optionally be used for on node GPU support instead cuBLAS/cuSOLVER, and performs slightly better.
1084*29e0a805SPieter Ghysels   ParMETIS and PTScotch can be used for parallel fill-reducing ordering.
1085*29e0a805SPieter Ghysels   ZFP is used for floating point compression of the sparse factors (LOSSY or LOSSLESS compression).
1086*29e0a805SPieter Ghysels   ButterflyPACK is used for HODLR (Hierarchically Off-Diagonal Low Rank) and HODBF (Hierarchically Off-Diagonal Butterfly) compression of the sparse factors.
1087*29e0a805SPieter Ghysels 
1088*29e0a805SPieter 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 
1092*29e0a805SPieter 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:
1097*29e0a805SPieter Ghysels + -mat_strumpack_verbose                      - Enable verbose output
1098*29e0a805SPieter 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
1099*29e0a805SPieter Ghysels . -mat_strumpack_compression_rel_tol          - Relative compression tolerance, when using `-pctype ilu`
1100*29e0a805SPieter Ghysels . -mat_strumpack_compression_abs_tol          - Absolute compression tolerance, when using `-pctype ilu`
1101*29e0a805SPieter Ghysels . -mat_strumpack_compression_min_sep_size     - Minimum size of separator for rank-structured compression, when using `-pctype ilu`
1102*29e0a805SPieter Ghysels . -mat_strumpack_compression_leaf_size        - Size of diagonal blocks in rank-structured approximation, when using `-pctype ilu`
1103*29e0a805SPieter Ghysels . -mat_strumpack_compression_lossy_precision  - Precision when using lossy compression [1-64], when using `-pctype ilu`, compression LOSSY (requires ZFP support)
1104*29e0a805SPieter 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)
1105*29e0a805SPieter Ghysels . -mat_strumpack_gpu                          - Enable GPU acceleration in numerical factorization (not supported for all compression types)
1106*29e0a805SPieter Ghysels . -mat_strumpack_colperm <TRUE>               - Permute matrix to make diagonal nonzeros
1107*29e0a805SPieter Ghysels . -mat_strumpack_reordering <METIS>           - Sparsity reducing matrix reordering (choose one of) NATURAL METIS PARMETIS SCOTCH PTSCOTCH RCM GEOMETRIC AMD MMD AND MLF SPECTRAL
1108*29e0a805SPieter Ghysels . -mat_strumpack_geometric_xyz <1,1,1>        - Mesh x,y,z dimensions, for use with GEOMETRIC ordering
1109*29e0a805SPieter Ghysels . -mat_strumpack_geometric_components <1>     - Number of components per mesh point, for geometric nested dissection ordering
1110*29e0a805SPieter Ghysels . -mat_strumpack_geometric_width <1>          - Width of the separator of the mesh, for geometric nested dissection ordering
1111*29e0a805SPieter 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 
1115*29e0a805SPieter Ghysels  HODLR, BLR_HODBF and ZFP_BLR_HODLR compression require STRUMPACK to be configured with ButterflyPACK support (`--download-butterflypack`).
1116*29e0a805SPieter Ghysels 
1117*29e0a805SPieter Ghysels  LOSSY, LOSSLESS and ZFP_BLR_HODLR compression require STRUMPACK to be configured with ZFP support (`--download-zfp`).
1118*29e0a805SPieter Ghysels 
11191cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `PCLU`, `PCILU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`, `PCFactorSetMatSolverType()`, `MatSolverType`,
1120*29e0a805SPieter 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;
1126*29e0a805SPieter Ghysels   PetscBool verb, flg, set;
1127*29e0a805SPieter Ghysels   PetscReal ctol;
1128*29e0a805SPieter Ghysels   PetscInt  min_sep_size, leaf_size, nxyz[3], nrdims, nc, w;
1129*29e0a805SPieter Ghysels #if defined(STRUMPACK_USE_ZFP)
1130*29e0a805SPieter Ghysels   PetscInt lossy_prec;
1131*29e0a805SPieter Ghysels #endif
1132*29e0a805SPieter Ghysels #if defined(STRUMPACK_USE_BPACK)
1133*29e0a805SPieter Ghysels   PetscInt bfly_lvls;
1134*29e0a805SPieter Ghysels #endif
1135*29e0a805SPieter Ghysels #if defined(STRUMPACK_USE_SLATE_SCALAPACK)
1136*29e0a805SPieter Ghysels   PetscMPIInt mpithreads;
1137*29e0a805SPieter Ghysels #endif
113834f43fa5SPieter Ghysels   STRUMPACK_SparseSolver       *S;
113934f43fa5SPieter Ghysels   STRUMPACK_INTERFACE           iface;
1140*29e0a805SPieter Ghysels   STRUMPACK_REORDERING_STRATEGY ndcurrent, ndvalue;
1141*29e0a805SPieter 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];
1147*29e0a805SPieter Ghysels   const char *const         STRUMPACKNDTypes[] = {"NATURAL", "METIS", "PARMETIS", "SCOTCH", "PTSCOTCH", "RCM", "GEOMETRIC", "AMD", "MMD", "AND", "MLF", "SPECTRAL", "STRUMPACKNDTypes", "", 0};
1148*29e0a805SPieter Ghysels   const char *const         CompTypes[]        = {"NONE", "HSS", "BLR", "HODLR", "BLR_HODLR", "ZFP_BLR_HODLR", "LOSSLESS", "LOSSY", "CompTypes", "", 0};
11497d6ea485SPieter Ghysels 
11507d6ea485SPieter Ghysels   PetscFunctionBegin;
1151*29e0a805SPieter Ghysels #if defined(STRUMPACK_USE_SLATE_SCALAPACK)
1152*29e0a805SPieter Ghysels   PetscCallMPI(MPI_Query_thread(&mpithreads));
1153*29e0a805SPieter Ghysels   PetscCheck(mpithreads == MPI_THREAD_MULTIPLE, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP_SYS, "SLATE requires MPI_THREAD_MULTIPLE");
1154*29e0a805SPieter 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));
1160*29e0a805SPieter Ghysels   PetscCall(MatSeqAIJSetPreallocation(B, 0, NULL));
1161*29e0a805SPieter 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));
1173*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetReordering_C", MatSTRUMPACKGetReordering_STRUMPACK));
11749566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetColPerm_C", MatSTRUMPACKSetColPerm_STRUMPACK));
1175*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetColPerm_C", MatSTRUMPACKGetColPerm_STRUMPACK));
1176*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetGeometricNxyz_C", MatSTRUMPACKSetGeometricNxyz_STRUMPACK));
1177*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetGeometricComponents_C", MatSTRUMPACKSetGeometricComponents_STRUMPACK));
1178*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetGeometricWidth_C", MatSTRUMPACKSetGeometricWidth_STRUMPACK));
1179*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetGPU_C", MatSTRUMPACKSetGPU_STRUMPACK));
1180*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetGPU_C", MatSTRUMPACKGetGPU_STRUMPACK));
1181*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetCompression_C", MatSTRUMPACKSetCompression_STRUMPACK));
1182*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetCompression_C", MatSTRUMPACKGetCompression_STRUMPACK));
1183*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetCompRelTol_C", MatSTRUMPACKSetCompRelTol_STRUMPACK));
1184*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetCompRelTol_C", MatSTRUMPACKGetCompRelTol_STRUMPACK));
1185*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetCompAbsTol_C", MatSTRUMPACKSetCompAbsTol_STRUMPACK));
1186*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetCompAbsTol_C", MatSTRUMPACKGetCompAbsTol_STRUMPACK));
1187*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetCompLeafSize_C", MatSTRUMPACKSetCompLeafSize_STRUMPACK));
1188*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetCompLeafSize_C", MatSTRUMPACKGetCompLeafSize_STRUMPACK));
1189*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetCompMinSepSize_C", MatSTRUMPACKSetCompMinSepSize_STRUMPACK));
1190*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetCompMinSepSize_C", MatSTRUMPACKGetCompMinSepSize_STRUMPACK));
1191*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetCompLossyPrecision_C", MatSTRUMPACKSetCompLossyPrecision_STRUMPACK));
1192*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetCompLossyPrecision_C", MatSTRUMPACKGetCompLossyPrecision_STRUMPACK));
1193*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetCompButterflyLevels_C", MatSTRUMPACKSetCompButterflyLevels_STRUMPACK));
1194*29e0a805SPieter Ghysels   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKGetCompButterflyLevels_C", MatSTRUMPACKGetCompButterflyLevels_STRUMPACK));
1195575704cbSPieter Ghysels   B->factortype = ftype;
1196*29e0a805SPieter Ghysels 
1197*29e0a805SPieter Ghysels   /* set solvertype */
1198*29e0a805SPieter Ghysels   PetscCall(PetscFree(B->solvertype));
1199*29e0a805SPieter Ghysels   PetscCall(PetscStrallocpy(MATSOLVERSTRUMPACK, &B->solvertype));
1200*29e0a805SPieter 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 
1213*29e0a805SPieter Ghysels   /* By default, no compression is done. Compression is enabled when the user enables it with        */
1214*29e0a805SPieter Ghysels   /*  -mat_strumpack_compression with anything else than NONE, or when selecting ilu                 */
1215*29e0a805SPieter Ghysels   /* preconditioning, in which case we default to STRUMPACK_BLR compression.                         */
1216*29e0a805SPieter Ghysels   /* When compression is enabled, the STRUMPACK solver becomes an incomplete                         */
121788382b05SPieter Ghysels   /* (or approximate) LU factorization.                                                              */
1218*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_compression", compcurrent = STRUMPACK_compression(*S));
1219*29e0a805SPieter Ghysels   PetscCall(PetscOptionsEnum("-mat_strumpack_compression", "Rank-structured compression type", "None", CompTypes, (PetscEnum)compcurrent, (PetscEnum *)&compvalue, &set));
1220*29e0a805SPieter Ghysels   if (set) {
1221*29e0a805SPieter Ghysels     PetscCall(MatSTRUMPACKSetCompression(B, (MatSTRUMPACKCompressionType)compvalue));
1222*29e0a805SPieter Ghysels   } else {
1223*29e0a805SPieter Ghysels     if (ftype == MAT_FACTOR_ILU) { PetscStackCallExternalVoid("STRUMPACK_set_compression", STRUMPACK_set_compression(*S, STRUMPACK_BLR)); }
122488382b05SPieter Ghysels   }
122555c022e5SPieter Ghysels 
1226*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_compression_rel_tol", ctol = (PetscReal)STRUMPACK_compression_rel_tol(*S));
1227*29e0a805SPieter Ghysels   PetscCall(PetscOptionsReal("-mat_strumpack_compression_rel_tol", "Relative compression tolerance", "None", ctol, &ctol, &set));
1228*29e0a805SPieter Ghysels   if (set) PetscStackCallExternalVoid("STRUMPACK_set_compression_rel_tol", STRUMPACK_set_compression_rel_tol(*S, (double)ctol));
1229*29e0a805SPieter Ghysels 
1230*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_compression_abs_tol", ctol = (PetscReal)STRUMPACK_compression_abs_tol(*S));
1231*29e0a805SPieter Ghysels   PetscCall(PetscOptionsReal("-mat_strumpack_compression_abs_tol", "Absolute compression tolerance", "None", ctol, &ctol, &set));
1232*29e0a805SPieter Ghysels   if (set) PetscStackCallExternalVoid("STRUMPACK_set_compression_abs_tol", STRUMPACK_set_compression_abs_tol(*S, (double)ctol));
1233*29e0a805SPieter Ghysels 
1234*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_compression_min_sep_size", min_sep_size = (PetscInt)STRUMPACK_compression_min_sep_size(*S));
1235*29e0a805SPieter Ghysels   PetscCall(PetscOptionsInt("-mat_strumpack_compression_min_sep_size", "Minimum size of separator for compression", "None", min_sep_size, &min_sep_size, &set));
1236*29e0a805SPieter Ghysels   if (set) PetscStackCallExternalVoid("STRUMPACK_set_compression_min_sep_size", STRUMPACK_set_compression_min_sep_size(*S, (int)min_sep_size));
1237*29e0a805SPieter Ghysels 
1238*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_compression_leaf_size", leaf_size = (PetscInt)STRUMPACK_compression_leaf_size(*S));
1239*29e0a805SPieter Ghysels   PetscCall(PetscOptionsInt("-mat_strumpack_compression_leaf_size", "Size of diagonal blocks in rank-structured approximation", "None", leaf_size, &leaf_size, &set));
1240*29e0a805SPieter Ghysels   if (set) PetscStackCallExternalVoid("STRUMPACK_set_compression_leaf_size", STRUMPACK_set_compression_leaf_size(*S, (int)leaf_size));
1241*29e0a805SPieter Ghysels 
1242*29e0a805SPieter Ghysels #if defined(STRUMPACK_USE_ZFP)
1243*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_compression_lossy_precision", lossy_prec = (PetscInt)STRUMPACK_compression_lossy_precision(*S));
1244*29e0a805SPieter Ghysels   PetscCall(PetscOptionsInt("-mat_strumpack_compression_lossy_precision", "Number of bitplanes to use in lossy compression", "None", lossy_prec, &lossy_prec, &set));
1245*29e0a805SPieter Ghysels   if (set) PetscStackCallExternalVoid("STRUMPACK_set_compression_lossy_precision", STRUMPACK_set_compression_lossy_precision(*S, (int)lossy_prec));
1246*29e0a805SPieter Ghysels #endif
1247*29e0a805SPieter Ghysels 
1248*29e0a805SPieter Ghysels #if defined(STRUMPACK_USE_BPACK)
1249*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_compression_butterfly_levels", bfly_lvls = (PetscInt)STRUMPACK_compression_butterfly_levels(*S));
1250*29e0a805SPieter 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));
1251*29e0a805SPieter Ghysels   if (set) PetscStackCallExternalVoid("STRUMPACK_set_compression_butterfly_levels", STRUMPACK_set_compression_butterfly_levels(*S, (int)bfly_lvls));
1252*29e0a805SPieter Ghysels #endif
1253*29e0a805SPieter Ghysels 
1254*29e0a805SPieter Ghysels #if defined(STRUMPACK_USE_CUDA) || defined(STRUMPACK_USE_HIP) || defined(STRUMPACK_USE_SYCL)
1255*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_use_gpu", flg = (STRUMPACK_use_gpu(*S) == 0) ? PETSC_FALSE : PETSC_TRUE);
1256*29e0a805SPieter Ghysels   PetscCall(PetscOptionsBool("-mat_strumpack_gpu", "Enable GPU acceleration (not supported for all compression types)", "None", flg, &flg, &set));
1257*29e0a805SPieter Ghysels   if (set) MatSTRUMPACKSetGPU(B, flg);
1258*29e0a805SPieter Ghysels #endif
1259*29e0a805SPieter Ghysels 
1260*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_matching", flg = (STRUMPACK_matching(*S) == STRUMPACK_MATCHING_NONE) ? PETSC_FALSE : PETSC_TRUE);
1261*29e0a805SPieter Ghysels   PetscCall(PetscOptionsBool("-mat_strumpack_colperm", "Find a col perm to get nonzero diagonal", "None", flg, &flg, &set));
1262*29e0a805SPieter Ghysels   if (set) PetscStackCallExternalVoid("STRUMPACK_set_matching", STRUMPACK_set_matching(*S, flg ? STRUMPACK_MATCHING_MAX_DIAGONAL_PRODUCT_SCALING : STRUMPACK_MATCHING_NONE));
1263*29e0a805SPieter Ghysels 
1264*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_reordering_method", ndcurrent = STRUMPACK_reordering_method(*S));
1265*29e0a805SPieter Ghysels   PetscCall(PetscOptionsEnum("-mat_strumpack_reordering", "Sparsity reducing matrix reordering", "None", STRUMPACKNDTypes, (PetscEnum)ndcurrent, (PetscEnum *)&ndvalue, &set));
1266*29e0a805SPieter Ghysels   if (set) PetscStackCallExternalVoid("STRUMPACK_set_reordering_method", STRUMPACK_set_reordering_method(*S, ndvalue));
1267*29e0a805SPieter Ghysels 
1268*29e0a805SPieter Ghysels   /* geometric ordering, for a regular 1D/2D/3D mesh in the natural ordering, */
1269*29e0a805SPieter Ghysels   /* with nc DOF's per gridpoint, and possibly a wider stencil                */
1270*29e0a805SPieter Ghysels   nrdims  = 3;
1271*29e0a805SPieter Ghysels   nxyz[0] = nxyz[1] = nxyz[2] = 1;
1272*29e0a805SPieter Ghysels   PetscCall(PetscOptionsIntArray("-mat_strumpack_geometric_xyz", "Mesh sizes nx,ny,nz (Use 1 for default)", "", nxyz, &nrdims, &set));
1273*29e0a805SPieter Ghysels   if (set) {
1274*29e0a805SPieter 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.");
1275*29e0a805SPieter Ghysels     PetscCall(MatSTRUMPACKSetGeometricNxyz(B, (int)nxyz[0], (int)nxyz[1], (int)nxyz[2]));
1276*29e0a805SPieter Ghysels   }
1277*29e0a805SPieter Ghysels   PetscCall(PetscOptionsInt("-mat_strumpack_geometric_components", "Number of components per mesh point, for geometric nested dissection ordering", "None", 1, &nc, &set));
1278*29e0a805SPieter Ghysels   if (set) PetscStackCallExternalVoid("STRUMPACK_set_components", STRUMPACK_set_components(*S, (int)nc));
1279*29e0a805SPieter 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));
1280*29e0a805SPieter Ghysels   if (set) PetscStackCallExternalVoid("STRUMPACK_set_separator_width", STRUMPACK_set_separator_width(*S, (int)w));
1281*29e0a805SPieter Ghysels 
1282*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_use_METIS_NodeNDP", flg = (STRUMPACK_use_METIS_NodeNDP(*S) == 0) ? PETSC_FALSE : PETSC_TRUE);
1283*29e0a805SPieter Ghysels   PetscCall(PetscOptionsBool("-mat_strumpack_metis_nodeNDP", "Use METIS_NodeNDP instead of METIS_NodeND, for a more balanced tree", "None", flg, &flg, &set));
1284*29e0a805SPieter Ghysels   if (set) {
1285*29e0a805SPieter Ghysels     if (flg) {
1286*29e0a805SPieter Ghysels       PetscStackCallExternalVoid("STRUMPACK_enable_METIS_NodeNDP", STRUMPACK_enable_METIS_NodeNDP(*S));
1287*29e0a805SPieter Ghysels     } else {
1288*29e0a805SPieter Ghysels       PetscStackCallExternalVoid("STRUMPACK_disable_METIS_NodeNDP", STRUMPACK_disable_METIS_NodeNDP(*S));
1289*29e0a805SPieter Ghysels     }
1290*29e0a805SPieter Ghysels   }
1291*29e0a805SPieter Ghysels 
1292*29e0a805SPieter Ghysels   /* Disable the outer iterative solver from STRUMPACK.                                       */
1293*29e0a805SPieter Ghysels   /* When STRUMPACK is used as a direct solver, it will by default do iterative refinement.   */
1294*29e0a805SPieter Ghysels   /* When STRUMPACK is used as an approximate factorization preconditioner (by enabling       */
1295*29e0a805SPieter Ghysels   /* low-rank compression), it will use it's own preconditioned GMRES. Here we can disable    */
1296*29e0a805SPieter Ghysels   /* the outer iterative solver, as PETSc uses STRUMPACK from within a KSP.                   */
1297*29e0a805SPieter Ghysels   PetscStackCallExternalVoid("STRUMPACK_set_Krylov_solver", STRUMPACK_set_Krylov_solver(*S, STRUMPACK_DIRECT));
1298*29e0a805SPieter Ghysels 
1299*29e0a805SPieter 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