xref: /petsc/src/mat/impls/aij/mpi/strumpack/strumpack.c (revision 2ef1f0ff6e3530e8731eb06ad663081f5844f49f)
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));
249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetColPerm_C", NULL));
259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetHSSRelTol_C", NULL));
269566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetHSSAbsTol_C", NULL));
279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetHSSMaxRank_C", NULL));
289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetHSSLeafSize_C", NULL));
299566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetHSSMinSepSize_C", NULL));
30575704cbSPieter Ghysels 
313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
32575704cbSPieter Ghysels }
33575704cbSPieter Ghysels 
34d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSTRUMPACKSetReordering_STRUMPACK(Mat F, MatSTRUMPACKReordering reordering)
35d71ae5a4SJacob Faibussowitsch {
3635e8bcc9SJunchao Zhang   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
37575704cbSPieter Ghysels 
38575704cbSPieter Ghysels   PetscFunctionBegin;
39e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_reordering_method", STRUMPACK_set_reordering_method(*S, (STRUMPACK_REORDERING_STRATEGY)reordering));
403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4134f43fa5SPieter Ghysels }
42e5a36eccSStefano Zampini 
43542aee0fSPieter Ghysels /*@
4434f43fa5SPieter Ghysels   MatSTRUMPACKSetReordering - Set STRUMPACK fill-reducing reordering
4534f43fa5SPieter Ghysels 
4634f43fa5SPieter Ghysels    Input Parameters:
47*2ef1f0ffSBarry Smith +  F - the factored matrix obtained by calling `MatGetFactor()`
48*2ef1f0ffSBarry Smith -  reordering - the code to be used to find the fill-reducing reordering, see `MatSTRUMPACKReordering`
4934f43fa5SPieter Ghysels 
5011a5261eSBarry Smith   Options Database Key:
51*2ef1f0ffSBarry Smith .   -mat_strumpack_reordering <METIS>  - Sparsity reducing matrix reordering, see `MatSTRUMPACKReordering`
5234f43fa5SPieter Ghysels 
5334f43fa5SPieter Ghysels    Level: beginner
5434f43fa5SPieter Ghysels 
5534f43fa5SPieter Ghysels    References:
56606c0280SSatish Balay .  * - STRUMPACK manual
5734f43fa5SPieter Ghysels 
58*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MatSTRUMPACKReordering`, `MatGetFactor()`, `MatSTRUMPACKSetColPerm()`, `MatSTRUMPACKSetHSSRelTol()`,
59*2ef1f0ffSBarry Smith           `MatSTRUMPACKSetHSSAbsTol()`, `MatSTRUMPACKSetHSSMaxRank()`, `MatSTRUMPACKSetHSSLeafSize()`, `MatSTRUMPACKSetHSSMinSepSize()`
60542aee0fSPieter Ghysels @*/
61d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSTRUMPACKSetReordering(Mat F, MatSTRUMPACKReordering reordering)
62d71ae5a4SJacob Faibussowitsch {
6334f43fa5SPieter Ghysels   PetscFunctionBegin;
6434f43fa5SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
6534f43fa5SPieter Ghysels   PetscValidLogicalCollectiveEnum(F, reordering, 2);
66cac4c232SBarry Smith   PetscTryMethod(F, "MatSTRUMPACKSetReordering_C", (Mat, MatSTRUMPACKReordering), (F, reordering));
673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
68575704cbSPieter Ghysels }
69575704cbSPieter Ghysels 
70d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSTRUMPACKSetColPerm_STRUMPACK(Mat F, PetscBool cperm)
71d71ae5a4SJacob Faibussowitsch {
7235e8bcc9SJunchao Zhang   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
7334f43fa5SPieter Ghysels 
7434f43fa5SPieter Ghysels   PetscFunctionBegin;
75e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_set_mc64job", STRUMPACK_set_mc64job(*S, cperm ? 5 : 0));
763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7734f43fa5SPieter Ghysels }
78e5a36eccSStefano Zampini 
79575704cbSPieter Ghysels /*@
8034f43fa5SPieter Ghysels   MatSTRUMPACKSetColPerm - Set whether STRUMPACK should try to permute the columns of the matrix in order to get a nonzero diagonal
81147403d9SBarry Smith 
82c3339decSBarry Smith    Logically Collective
83575704cbSPieter Ghysels 
84575704cbSPieter Ghysels    Input Parameters:
85*2ef1f0ffSBarry Smith +  F - the factored matrix obtained by calling `MatGetFactor()`
8611a5261eSBarry Smith -  cperm - `PETSC_TRUE` to permute (internally) the columns of the matrix
87575704cbSPieter Ghysels 
8811a5261eSBarry Smith   Options Database Key:
89147403d9SBarry Smith .   -mat_strumpack_colperm <cperm> - true to use the permutation
90575704cbSPieter Ghysels 
91575704cbSPieter Ghysels    Level: beginner
92575704cbSPieter Ghysels 
93575704cbSPieter Ghysels    References:
94606c0280SSatish Balay .  * - STRUMPACK manual
95575704cbSPieter Ghysels 
96*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `MatSTRUMPACKSetReordering()`, `Mat`, `MatGetFactor()`, `MatSTRUMPACKSetHSSRelTol()`, `MatSTRUMPACKSetHSSAbsTol()`,
97*2ef1f0ffSBarry Smith           `MatSTRUMPACKSetHSSMaxRank()`, `MatSTRUMPACKSetHSSLeafSize()`, `MatSTRUMPACKSetHSSMinSepSize()`
98575704cbSPieter Ghysels @*/
99d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSTRUMPACKSetColPerm(Mat F, PetscBool cperm)
100d71ae5a4SJacob Faibussowitsch {
101575704cbSPieter Ghysels   PetscFunctionBegin;
102575704cbSPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
10334f43fa5SPieter Ghysels   PetscValidLogicalCollectiveBool(F, cperm, 2);
104cac4c232SBarry Smith   PetscTryMethod(F, "MatSTRUMPACKSetColPerm_C", (Mat, PetscBool), (F, cperm));
1053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
106575704cbSPieter Ghysels }
107575704cbSPieter Ghysels 
108d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSTRUMPACKSetHSSRelTol_STRUMPACK(Mat F, PetscReal rtol)
109d71ae5a4SJacob Faibussowitsch {
11035e8bcc9SJunchao Zhang   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
11134f43fa5SPieter Ghysels 
11234f43fa5SPieter Ghysels   PetscFunctionBegin;
113e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_set_HSS_rel_tol", STRUMPACK_set_HSS_rel_tol(*S, rtol));
1143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11534f43fa5SPieter Ghysels }
116e5a36eccSStefano Zampini 
11734f43fa5SPieter Ghysels /*@
11834f43fa5SPieter Ghysels   MatSTRUMPACKSetHSSRelTol - Set STRUMPACK relative tolerance for HSS compression
119147403d9SBarry Smith 
120c3339decSBarry Smith   Logically Collective
12134f43fa5SPieter Ghysels 
12234f43fa5SPieter Ghysels    Input Parameters:
123*2ef1f0ffSBarry Smith +  F - the factored matrix obtained by calling `MatGetFactor()`
12434f43fa5SPieter Ghysels -  rtol - relative compression tolerance
12534f43fa5SPieter Ghysels 
12611a5261eSBarry Smith   Options Database Key:
12734f43fa5SPieter Ghysels .   -mat_strumpack_hss_rel_tol <1e-2>         - Relative compression tolerance (None)
12834f43fa5SPieter Ghysels 
12934f43fa5SPieter Ghysels    Level: beginner
13034f43fa5SPieter Ghysels 
13134f43fa5SPieter Ghysels    References:
132606c0280SSatish Balay .  * - STRUMPACK manual
13334f43fa5SPieter Ghysels 
134*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MatGetFactor()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()`, `MatSTRUMPACKSetHSSAbsTol()`,
135*2ef1f0ffSBarry Smith           `MatSTRUMPACKSetHSSMaxRank()`, `MatSTRUMPACKSetHSSLeafSize()`, `MatSTRUMPACKSetHSSMinSepSize()`
13634f43fa5SPieter Ghysels @*/
137d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSTRUMPACKSetHSSRelTol(Mat F, PetscReal rtol)
138d71ae5a4SJacob Faibussowitsch {
13934f43fa5SPieter Ghysels   PetscFunctionBegin;
14034f43fa5SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
14134f43fa5SPieter Ghysels   PetscValidLogicalCollectiveReal(F, rtol, 2);
142cac4c232SBarry Smith   PetscTryMethod(F, "MatSTRUMPACKSetHSSRelTol_C", (Mat, PetscReal), (F, rtol));
1433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14434f43fa5SPieter Ghysels }
14534f43fa5SPieter Ghysels 
146d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSTRUMPACKSetHSSAbsTol_STRUMPACK(Mat F, PetscReal atol)
147d71ae5a4SJacob Faibussowitsch {
14835e8bcc9SJunchao Zhang   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
14934f43fa5SPieter Ghysels 
15034f43fa5SPieter Ghysels   PetscFunctionBegin;
151e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_set_HSS_abs_tol", STRUMPACK_set_HSS_abs_tol(*S, atol));
1523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15334f43fa5SPieter Ghysels }
154e5a36eccSStefano Zampini 
15534f43fa5SPieter Ghysels /*@
15634f43fa5SPieter Ghysels   MatSTRUMPACKSetHSSAbsTol - Set STRUMPACK absolute tolerance for HSS compression
157147403d9SBarry Smith 
158c3339decSBarry Smith    Logically Collective
15934f43fa5SPieter Ghysels 
16034f43fa5SPieter Ghysels    Input Parameters:
161*2ef1f0ffSBarry Smith +  F - the factored matrix obtained by calling `MatGetFactor()`
16234f43fa5SPieter Ghysels -  atol - absolute compression tolerance
16334f43fa5SPieter Ghysels 
16411a5261eSBarry Smith   Options Database Key:
16534f43fa5SPieter Ghysels .   -mat_strumpack_hss_abs_tol <1e-8>         - Absolute compression tolerance (None)
16634f43fa5SPieter Ghysels 
16734f43fa5SPieter Ghysels    Level: beginner
16834f43fa5SPieter Ghysels 
16934f43fa5SPieter Ghysels    References:
170606c0280SSatish Balay .  * - STRUMPACK manual
17134f43fa5SPieter Ghysels 
172*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MatGetFactor()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()`, `MatSTRUMPACKSetHSSRelTol()`,
173*2ef1f0ffSBarry Smith           `MatSTRUMPACKSetHSSMaxRank()`, `MatSTRUMPACKSetHSSLeafSize()`, `MatSTRUMPACKSetHSSMinSepSize()`
17434f43fa5SPieter Ghysels @*/
175d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSTRUMPACKSetHSSAbsTol(Mat F, PetscReal atol)
176d71ae5a4SJacob Faibussowitsch {
17734f43fa5SPieter Ghysels   PetscFunctionBegin;
17834f43fa5SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
17934f43fa5SPieter Ghysels   PetscValidLogicalCollectiveReal(F, atol, 2);
180cac4c232SBarry Smith   PetscTryMethod(F, "MatSTRUMPACKSetHSSAbsTol_C", (Mat, PetscReal), (F, atol));
1813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18234f43fa5SPieter Ghysels }
18334f43fa5SPieter Ghysels 
184d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSTRUMPACKSetHSSMaxRank_STRUMPACK(Mat F, PetscInt hssmaxrank)
185d71ae5a4SJacob Faibussowitsch {
18635e8bcc9SJunchao Zhang   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
18734f43fa5SPieter Ghysels 
18834f43fa5SPieter Ghysels   PetscFunctionBegin;
189e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_set_HSS_max_rank", STRUMPACK_set_HSS_max_rank(*S, hssmaxrank));
1903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19134f43fa5SPieter Ghysels }
192e5a36eccSStefano Zampini 
19334f43fa5SPieter Ghysels /*@
19434f43fa5SPieter Ghysels   MatSTRUMPACKSetHSSMaxRank - Set STRUMPACK maximum HSS rank
195147403d9SBarry Smith 
196c3339decSBarry Smith    Logically Collective
19734f43fa5SPieter Ghysels 
19834f43fa5SPieter Ghysels    Input Parameters:
199*2ef1f0ffSBarry Smith +  F - the factored matrix obtained by calling `MatGetFactor()`
20034f43fa5SPieter Ghysels -  hssmaxrank - maximum rank used in low-rank approximation
20134f43fa5SPieter Ghysels 
20211a5261eSBarry Smith   Options Database Key:
20334f43fa5SPieter Ghysels .   -mat_strumpack_max_rank - Maximum rank in HSS compression, when using pctype ilu (None)
20434f43fa5SPieter Ghysels 
20534f43fa5SPieter Ghysels    Level: beginner
20634f43fa5SPieter Ghysels 
20734f43fa5SPieter Ghysels    References:
208606c0280SSatish Balay .  * - STRUMPACK manual
20934f43fa5SPieter Ghysels 
210*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MatGetFactor()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()`, `MatSTRUMPACKSetHSSRelTol()`,
211*2ef1f0ffSBarry Smith           `MatSTRUMPACKSetHSSAbsTol()`, `MatSTRUMPACKSetHSSLeafSize()`, `MatSTRUMPACKSetHSSMinSepSize()`
21234f43fa5SPieter Ghysels @*/
213d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSTRUMPACKSetHSSMaxRank(Mat F, PetscInt hssmaxrank)
214d71ae5a4SJacob Faibussowitsch {
21534f43fa5SPieter Ghysels   PetscFunctionBegin;
21634f43fa5SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
21734f43fa5SPieter Ghysels   PetscValidLogicalCollectiveInt(F, hssmaxrank, 2);
218cac4c232SBarry Smith   PetscTryMethod(F, "MatSTRUMPACKSetHSSMaxRank_C", (Mat, PetscInt), (F, hssmaxrank));
2193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22034f43fa5SPieter Ghysels }
22134f43fa5SPieter Ghysels 
222d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSTRUMPACKSetHSSLeafSize_STRUMPACK(Mat F, PetscInt leaf_size)
223d71ae5a4SJacob Faibussowitsch {
22435e8bcc9SJunchao Zhang   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
225a36bf211SPieter Ghysels 
226a36bf211SPieter Ghysels   PetscFunctionBegin;
227e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_set_HSS_leaf_size", STRUMPACK_set_HSS_leaf_size(*S, leaf_size));
2283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
229a36bf211SPieter Ghysels }
230e5a36eccSStefano Zampini 
231a36bf211SPieter Ghysels /*@
232a36bf211SPieter Ghysels   MatSTRUMPACKSetHSSLeafSize - Set STRUMPACK HSS leaf size
233147403d9SBarry Smith 
234c3339decSBarry Smith    Logically Collective
235a36bf211SPieter Ghysels 
236a36bf211SPieter Ghysels    Input Parameters:
237*2ef1f0ffSBarry Smith +  F - the factored matrix obtained by calling `MatGetFactor()`
238a36bf211SPieter Ghysels -  leaf_size - Size of diagonal blocks in HSS approximation
239a36bf211SPieter Ghysels 
24011a5261eSBarry Smith   Options Database Key:
241a36bf211SPieter Ghysels .   -mat_strumpack_leaf_size - Size of diagonal blocks in HSS approximation, when using pctype ilu (None)
242a36bf211SPieter Ghysels 
243a36bf211SPieter Ghysels    Level: beginner
244a36bf211SPieter Ghysels 
245a36bf211SPieter Ghysels    References:
246606c0280SSatish Balay .  * - STRUMPACK manual
247a36bf211SPieter Ghysels 
248*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MatGetFactor()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()`, `MatSTRUMPACKSetHSSRelTol()`,
249*2ef1f0ffSBarry Smith           `MatSTRUMPACKSetHSSAbsTol()`, `MatSTRUMPACKSetHSSMaxRank()`, `MatSTRUMPACKSetHSSMinSepSize()`
250a36bf211SPieter Ghysels @*/
251d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSTRUMPACKSetHSSLeafSize(Mat F, PetscInt leaf_size)
252d71ae5a4SJacob Faibussowitsch {
253a36bf211SPieter Ghysels   PetscFunctionBegin;
254a36bf211SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
255a36bf211SPieter Ghysels   PetscValidLogicalCollectiveInt(F, leaf_size, 2);
256cac4c232SBarry Smith   PetscTryMethod(F, "MatSTRUMPACKSetHSSLeafSize_C", (Mat, PetscInt), (F, leaf_size));
2573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
258a36bf211SPieter Ghysels }
259a36bf211SPieter Ghysels 
260d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSTRUMPACKSetHSSMinSepSize_STRUMPACK(Mat F, PetscInt hssminsize)
261d71ae5a4SJacob Faibussowitsch {
26235e8bcc9SJunchao Zhang   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
26334f43fa5SPieter Ghysels 
26434f43fa5SPieter Ghysels   PetscFunctionBegin;
265e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_set_HSS_min_sep_size", STRUMPACK_set_HSS_min_sep_size(*S, hssminsize));
2663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
26734f43fa5SPieter Ghysels }
268e5a36eccSStefano Zampini 
269291d0ed5SPieter Ghysels /*@
270291d0ed5SPieter Ghysels   MatSTRUMPACKSetHSSMinSepSize - Set STRUMPACK minimum separator size for low-rank approximation
271147403d9SBarry Smith 
272c3339decSBarry Smith    Logically Collective
273291d0ed5SPieter Ghysels 
274291d0ed5SPieter Ghysels    Input Parameters:
275*2ef1f0ffSBarry Smith +  F - the factored matrix obtained by calling `MatGetFactor()`
276291d0ed5SPieter Ghysels -  hssminsize - minimum dense matrix size for low-rank approximation
277291d0ed5SPieter Ghysels 
27811a5261eSBarry Smith   Options Database Key:
279147403d9SBarry Smith .   -mat_strumpack_hss_min_sep_size <hssminsize> - set the minimum separator size
280291d0ed5SPieter Ghysels 
281291d0ed5SPieter Ghysels    Level: beginner
282291d0ed5SPieter Ghysels 
283291d0ed5SPieter Ghysels    References:
284606c0280SSatish Balay .  * - STRUMPACK manual
285291d0ed5SPieter Ghysels 
286*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MatGetFactor()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()`, `MatSTRUMPACKSetHSSRelTol()`,
287*2ef1f0ffSBarry Smith           `MatSTRUMPACKSetHSSAbsTol()`, `MatSTRUMPACKSetHSSMaxRank()`, `MatSTRUMPACKSetHSSLeafSize()`
288291d0ed5SPieter Ghysels @*/
289d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat F, PetscInt hssminsize)
290d71ae5a4SJacob Faibussowitsch {
291291d0ed5SPieter Ghysels   PetscFunctionBegin;
292291d0ed5SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
293291d0ed5SPieter Ghysels   PetscValidLogicalCollectiveInt(F, hssminsize, 2);
294cac4c232SBarry Smith   PetscTryMethod(F, "MatSTRUMPACKSetHSSMinSepSize_C", (Mat, PetscInt), (F, hssminsize));
2953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
296291d0ed5SPieter Ghysels }
297291d0ed5SPieter Ghysels 
298d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolve_STRUMPACK(Mat A, Vec b_mpi, Vec x)
299d71ae5a4SJacob Faibussowitsch {
30035e8bcc9SJunchao Zhang   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)A->data;
3017d6ea485SPieter Ghysels   STRUMPACK_RETURN_CODE   sp_err;
3027d6ea485SPieter Ghysels   const PetscScalar      *bptr;
3037d6ea485SPieter Ghysels   PetscScalar            *xptr;
3047d6ea485SPieter Ghysels 
3057d6ea485SPieter Ghysels   PetscFunctionBegin;
3069566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x, &xptr));
3079566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(b_mpi, &bptr));
3080d08a34dSPieter Ghysels 
309e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_solve", sp_err = STRUMPACK_solve(*S, (PetscScalar *)bptr, xptr, 0));
3100d08a34dSPieter Ghysels   switch (sp_err) {
311d71ae5a4SJacob Faibussowitsch   case STRUMPACK_SUCCESS:
312d71ae5a4SJacob Faibussowitsch     break;
3139371c9d4SSatish Balay   case STRUMPACK_MATRIX_NOT_SET: {
3149371c9d4SSatish Balay     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix was not set");
3159371c9d4SSatish Balay     break;
3169371c9d4SSatish Balay   }
3179371c9d4SSatish Balay   case STRUMPACK_REORDERING_ERROR: {
3189371c9d4SSatish Balay     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix reordering failed");
3199371c9d4SSatish Balay     break;
3209371c9d4SSatish Balay   }
321d71ae5a4SJacob Faibussowitsch   default:
322d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: solve failed");
3237d6ea485SPieter Ghysels   }
3249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x, &xptr));
3259566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(b_mpi, &bptr));
3263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3277d6ea485SPieter Ghysels }
3287d6ea485SPieter Ghysels 
329d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMatSolve_STRUMPACK(Mat A, Mat B_mpi, Mat X)
330d71ae5a4SJacob Faibussowitsch {
3317d6ea485SPieter Ghysels   PetscBool flg;
3327d6ea485SPieter Ghysels 
3337d6ea485SPieter Ghysels   PetscFunctionBegin;
3349566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)B_mpi, &flg, MATSEQDENSE, MATMPIDENSE, NULL));
3355f80ce2aSJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Matrix B must be MATDENSE matrix");
3369566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)X, &flg, MATSEQDENSE, MATMPIDENSE, NULL));
3375f80ce2aSJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Matrix X must be MATDENSE matrix");
3387d6ea485SPieter Ghysels   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "MatMatSolve_STRUMPACK() is not implemented yet");
3393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3407d6ea485SPieter Ghysels }
3417d6ea485SPieter Ghysels 
342d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_Info_STRUMPACK(Mat A, PetscViewer viewer)
343d71ae5a4SJacob Faibussowitsch {
344ad0c5e61SPieter Ghysels   PetscFunctionBegin;
345ad0c5e61SPieter Ghysels   /* check if matrix is strumpack type */
3463ba16761SJacob Faibussowitsch   if (A->ops->solve != MatSolve_STRUMPACK) PetscFunctionReturn(PETSC_SUCCESS);
3479566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "STRUMPACK sparse solver!\n"));
3483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
349ad0c5e61SPieter Ghysels }
350ad0c5e61SPieter Ghysels 
351d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_STRUMPACK(Mat A, PetscViewer viewer)
352d71ae5a4SJacob Faibussowitsch {
353ad0c5e61SPieter Ghysels   PetscBool         iascii;
354ad0c5e61SPieter Ghysels   PetscViewerFormat format;
355ad0c5e61SPieter Ghysels 
356ad0c5e61SPieter Ghysels   PetscFunctionBegin;
3579566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
358ad0c5e61SPieter Ghysels   if (iascii) {
3599566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
3609566063dSJacob Faibussowitsch     if (format == PETSC_VIEWER_ASCII_INFO) PetscCall(MatView_Info_STRUMPACK(A, viewer));
361ad0c5e61SPieter Ghysels   }
3623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
363ad0c5e61SPieter Ghysels }
3647d6ea485SPieter Ghysels 
365d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLUFactorNumeric_STRUMPACK(Mat F, Mat A, const MatFactorInfo *info)
366d71ae5a4SJacob Faibussowitsch {
36735e8bcc9SJunchao Zhang   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->data;
3687d6ea485SPieter Ghysels   STRUMPACK_RETURN_CODE   sp_err;
3697d6ea485SPieter Ghysels   Mat_SeqAIJ             *A_d, *A_o;
3707d6ea485SPieter Ghysels   Mat_MPIAIJ             *mat;
3710d08a34dSPieter Ghysels   PetscInt                M = A->rmap->N, m = A->rmap->n;
3727d6ea485SPieter Ghysels   PetscBool               flg;
37335e8bcc9SJunchao Zhang   const PetscScalar      *A_d_a, *A_o_a;
3747d6ea485SPieter Ghysels 
3757d6ea485SPieter Ghysels   PetscFunctionBegin;
37635e8bcc9SJunchao Zhang   PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATMPIAIJ, &flg));
37735e8bcc9SJunchao Zhang   if (flg) { /* A might be MATMPIAIJCUSPARSE etc */
3787d6ea485SPieter Ghysels     mat = (Mat_MPIAIJ *)A->data;
37935e8bcc9SJunchao Zhang     PetscCall(MatSeqAIJGetArrayRead(mat->A, &A_d_a)); /* Make sure mat is sync'ed to host */
38035e8bcc9SJunchao Zhang     PetscCall(MatSeqAIJGetArrayRead(mat->B, &A_o_a));
3817d6ea485SPieter Ghysels     A_d = (Mat_SeqAIJ *)(mat->A)->data;
3827d6ea485SPieter Ghysels     A_o = (Mat_SeqAIJ *)(mat->B)->data;
38335e8bcc9SJunchao Zhang     PetscStackCallExternalVoid("STRUMPACK_set_MPIAIJ_matrix", STRUMPACK_set_MPIAIJ_matrix(*S, &m, A_d->i, A_d->j, A_d_a, A_o->i, A_o->j, A_o_a, mat->garray));
38435e8bcc9SJunchao Zhang     PetscCall(MatSeqAIJRestoreArrayRead(mat->A, &A_d_a));
38535e8bcc9SJunchao Zhang     PetscCall(MatSeqAIJRestoreArrayRead(mat->B, &A_o_a));
38635e8bcc9SJunchao Zhang   } else { /* A might be MATSEQAIJCUSPARSE etc */
38735e8bcc9SJunchao Zhang     PetscCall(MatSeqAIJGetArrayRead(A, &A_d_a));
3887d6ea485SPieter Ghysels     A_d = (Mat_SeqAIJ *)A->data;
38935e8bcc9SJunchao Zhang     PetscStackCallExternalVoid("STRUMPACK_set_csr_matrix", STRUMPACK_set_csr_matrix(*S, &M, A_d->i, A_d->j, A_d_a, 0));
39035e8bcc9SJunchao Zhang     PetscCall(MatSeqAIJRestoreArrayRead(A, &A_d_a));
3917d6ea485SPieter Ghysels   }
3927d6ea485SPieter Ghysels 
3937d6ea485SPieter Ghysels   /* Reorder and Factor the matrix. */
3947d6ea485SPieter Ghysels   /* TODO figure out how to avoid reorder if the matrix values changed, but the pattern remains the same. */
395e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_reorder", sp_err = STRUMPACK_reorder(*S));
396e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_factor", sp_err = STRUMPACK_factor(*S));
3970d08a34dSPieter Ghysels   switch (sp_err) {
398d71ae5a4SJacob Faibussowitsch   case STRUMPACK_SUCCESS:
399d71ae5a4SJacob Faibussowitsch     break;
4009371c9d4SSatish Balay   case STRUMPACK_MATRIX_NOT_SET: {
4019371c9d4SSatish Balay     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix was not set");
4029371c9d4SSatish Balay     break;
4039371c9d4SSatish Balay   }
4049371c9d4SSatish Balay   case STRUMPACK_REORDERING_ERROR: {
4059371c9d4SSatish Balay     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix reordering failed");
4069371c9d4SSatish Balay     break;
4079371c9d4SSatish Balay   }
408d71ae5a4SJacob Faibussowitsch   default:
409d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: factorization failed");
4107d6ea485SPieter Ghysels   }
411cb250fa3SPieter Ghysels   F->assembled    = PETSC_TRUE;
412cb250fa3SPieter Ghysels   F->preallocated = PETSC_TRUE;
4133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4147d6ea485SPieter Ghysels }
4157d6ea485SPieter Ghysels 
416d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLUFactorSymbolic_STRUMPACK(Mat F, Mat A, IS r, IS c, const MatFactorInfo *info)
417d71ae5a4SJacob Faibussowitsch {
41835e8bcc9SJunchao Zhang   STRUMPACK_SparseSolver       *S = (STRUMPACK_SparseSolver *)F->data;
41926cc229bSBarry Smith   PetscBool                     flg, set;
42026cc229bSBarry Smith   PetscReal                     ctol;
42126cc229bSBarry Smith   PetscInt                      hssminsize, max_rank, leaf_size;
42226cc229bSBarry Smith   STRUMPACK_REORDERING_STRATEGY ndcurrent, ndvalue;
42326cc229bSBarry Smith   STRUMPACK_KRYLOV_SOLVER       itcurrent, itsolver;
42426cc229bSBarry Smith   const char *const             STRUMPACKNDTypes[] = {"NATURAL", "METIS", "PARMETIS", "SCOTCH", "PTSCOTCH", "RCM", "STRUMPACKNDTypes", "", 0};
42526cc229bSBarry Smith   const char *const             SolverTypes[]      = {"AUTO", "NONE", "REFINE", "PREC_GMRES", "GMRES", "PREC_BICGSTAB", "BICGSTAB", "SolverTypes", "", 0};
42626cc229bSBarry Smith 
4277d6ea485SPieter Ghysels   PetscFunctionBegin;
42826cc229bSBarry Smith   /* Set options to F */
42926cc229bSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)F), ((PetscObject)F)->prefix, "STRUMPACK Options", "Mat");
430e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_HSS_rel_tol", ctol = (PetscReal)STRUMPACK_HSS_rel_tol(*S));
43126cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_strumpack_hss_rel_tol", "Relative compression tolerance", "None", ctol, &ctol, &set));
432e77caa6dSBarry Smith   if (set) PetscStackCallExternalVoid("STRUMPACK_set_HSS_rel_tol", STRUMPACK_set_HSS_rel_tol(*S, (double)ctol));
43326cc229bSBarry Smith 
434e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_HSS_abs_tol", ctol = (PetscReal)STRUMPACK_HSS_abs_tol(*S));
43526cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_strumpack_hss_abs_tol", "Absolute compression tolerance", "None", ctol, &ctol, &set));
436e77caa6dSBarry Smith   if (set) PetscStackCallExternalVoid("STRUMPACK_set_HSS_abs_tol", STRUMPACK_set_HSS_abs_tol(*S, (double)ctol));
43726cc229bSBarry Smith 
438e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_mc64job", flg = (STRUMPACK_mc64job(*S) == 0) ? PETSC_FALSE : PETSC_TRUE);
43926cc229bSBarry Smith   PetscCall(PetscOptionsBool("-mat_strumpack_colperm", "Find a col perm to get nonzero diagonal", "None", flg, &flg, &set));
440e77caa6dSBarry Smith   if (set) PetscStackCallExternalVoid("STRUMPACK_set_mc64job", STRUMPACK_set_mc64job(*S, flg ? 5 : 0));
44126cc229bSBarry Smith 
442e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_HSS_min_sep_size", hssminsize = (PetscInt)STRUMPACK_HSS_min_sep_size(*S));
44326cc229bSBarry Smith   PetscCall(PetscOptionsInt("-mat_strumpack_hss_min_sep_size", "Minimum size of separator for HSS compression", "None", hssminsize, &hssminsize, &set));
444e77caa6dSBarry Smith   if (set) PetscStackCallExternalVoid("STRUMPACK_set_HSS_min_sep_size", STRUMPACK_set_HSS_min_sep_size(*S, (int)hssminsize));
44526cc229bSBarry Smith 
446e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_HSS_max_rank", max_rank = (PetscInt)STRUMPACK_HSS_max_rank(*S));
44726cc229bSBarry Smith   PetscCall(PetscOptionsInt("-mat_strumpack_max_rank", "Maximum rank in HSS approximation", "None", max_rank, &max_rank, &set));
448e77caa6dSBarry Smith   if (set) PetscStackCallExternalVoid("STRUMPACK_set_HSS_max_rank", STRUMPACK_set_HSS_max_rank(*S, (int)max_rank));
44926cc229bSBarry Smith 
450e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_HSS_leaf_size", leaf_size = (PetscInt)STRUMPACK_HSS_leaf_size(*S));
45126cc229bSBarry Smith   PetscCall(PetscOptionsInt("-mat_strumpack_leaf_size", "Size of diagonal blocks in HSS approximation", "None", leaf_size, &leaf_size, &set));
452e77caa6dSBarry Smith   if (set) PetscStackCallExternalVoid("STRUMPACK_set_HSS_leaf_size", STRUMPACK_set_HSS_leaf_size(*S, (int)leaf_size));
45326cc229bSBarry Smith 
454e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_reordering_method", ndcurrent = STRUMPACK_reordering_method(*S));
45526cc229bSBarry Smith   PetscCall(PetscOptionsEnum("-mat_strumpack_reordering", "Sparsity reducing matrix reordering", "None", STRUMPACKNDTypes, (PetscEnum)ndcurrent, (PetscEnum *)&ndvalue, &set));
456e77caa6dSBarry Smith   if (set) PetscStackCallExternalVoid("STRUMPACK_set_reordering_method", STRUMPACK_set_reordering_method(*S, ndvalue));
45726cc229bSBarry Smith 
45826cc229bSBarry Smith   /* Disable the outer iterative solver from STRUMPACK.                                       */
45926cc229bSBarry Smith   /* When STRUMPACK is used as a direct solver, it will by default do iterative refinement.   */
46026cc229bSBarry Smith   /* When STRUMPACK is used as an approximate factorization preconditioner (by enabling       */
46126cc229bSBarry Smith   /* low-rank compression), it will use it's own preconditioned GMRES. Here we can disable    */
46226cc229bSBarry Smith   /* the outer iterative solver, as PETSc uses STRUMPACK from within a KSP.                   */
463e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_set_Krylov_solver", STRUMPACK_set_Krylov_solver(*S, STRUMPACK_DIRECT));
46426cc229bSBarry Smith 
465e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_Krylov_solver", itcurrent = STRUMPACK_Krylov_solver(*S));
46626cc229bSBarry Smith   PetscCall(PetscOptionsEnum("-mat_strumpack_iterative_solver", "Select iterative solver from STRUMPACK", "None", SolverTypes, (PetscEnum)itcurrent, (PetscEnum *)&itsolver, &set));
467e77caa6dSBarry Smith   if (set) PetscStackCallExternalVoid("STRUMPACK_set_Krylov_solver", STRUMPACK_set_Krylov_solver(*S, itsolver));
46826cc229bSBarry Smith   PetscOptionsEnd();
46926cc229bSBarry Smith 
4707d6ea485SPieter Ghysels   F->ops->lufactornumeric = MatLUFactorNumeric_STRUMPACK;
4717d6ea485SPieter Ghysels   F->ops->solve           = MatSolve_STRUMPACK;
4727d6ea485SPieter Ghysels   F->ops->matsolve        = MatMatSolve_STRUMPACK;
4733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4747d6ea485SPieter Ghysels }
4757d6ea485SPieter Ghysels 
476d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorGetSolverType_aij_strumpack(Mat A, MatSolverType *type)
477d71ae5a4SJacob Faibussowitsch {
4787d6ea485SPieter Ghysels   PetscFunctionBegin;
4797d6ea485SPieter Ghysels   *type = MATSOLVERSTRUMPACK;
4803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4817d6ea485SPieter Ghysels }
4827d6ea485SPieter Ghysels 
483575704cbSPieter Ghysels /*MC
484575704cbSPieter Ghysels   MATSOLVERSSTRUMPACK = "strumpack" - A solver package providing a direct sparse solver (PCLU)
485575704cbSPieter Ghysels   and a preconditioner (PCILU) using low-rank compression via the external package STRUMPACK.
486575704cbSPieter Ghysels 
487575704cbSPieter Ghysels   Consult the STRUMPACK-sparse manual for more info.
488575704cbSPieter Ghysels 
489*2ef1f0ffSBarry Smith   Use ` ./configure --download-strumpack` to have PETSc installed with STRUMPACK
490575704cbSPieter Ghysels 
491*2ef1f0ffSBarry Smith   Use `-pc_type lu` `-pc_factor_mat_solver_type strumpack` to use this as an exact (direct) solver.
492575704cbSPieter Ghysels 
493*2ef1f0ffSBarry Smith   Use `-pc_type ilu` `-pc_factor_mat_solver_type strumpack` to enable low-rank compression (i.e, use as a preconditioner).
494*2ef1f0ffSBarry Smith 
495*2ef1f0ffSBarry Smith   Works with `MATAIJ` matrices
496575704cbSPieter Ghysels 
497575704cbSPieter Ghysels   Options Database Keys:
49867b8a455SSatish Balay + -mat_strumpack_verbose                    - verbose info
49934f43fa5SPieter Ghysels . -mat_strumpack_hss_rel_tol <1e-2>         - Relative compression tolerance (None)
50034f43fa5SPieter Ghysels . -mat_strumpack_hss_abs_tol <1e-8>         - Absolute compression tolerance (None)
501575704cbSPieter Ghysels . -mat_strumpack_colperm <TRUE>             - Permute matrix to make diagonal nonzeros (None)
50234f43fa5SPieter Ghysels . -mat_strumpack_hss_min_sep_size <256>     - Minimum size of separator for HSS compression (None)
50334f43fa5SPieter Ghysels . -mat_strumpack_max_rank                   - Maximum rank in HSS compression, when using pctype ilu (None)
504a36bf211SPieter Ghysels . -mat_strumpack_leaf_size                  - Size of diagonal blocks in HSS approximation, when using pctype ilu (None)
505*2ef1f0ffSBarry Smith . -mat_strumpack_reordering <METIS>         - Sparsity reducing matrix reordering see `MatSTRUMPACKReordering`
506*2ef1f0ffSBarry Smith - -mat_strumpack_iterative_solver <DIRECT>  - Select iterative solver from STRUMPACK (choose one of) `AUTO`, `DIRECT`, `REFINE`, `PREC_GMRES`,
507*2ef1f0ffSBarry Smith                                               `GMRES`, `PREC_BICGSTAB`, `BICGSTAB`
508575704cbSPieter Ghysels 
509575704cbSPieter Ghysels  Level: beginner
510575704cbSPieter Ghysels 
511*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `PCLU`, `PCILU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`, `PCFactorSetMatSolverType()`, `MatSolverType`,
512*2ef1f0ffSBarry Smith           `MatGetFactor()`, `MatSTRUMPACKSetReordering()`, `MatSTRUMPACKSetColPerm()`, `MatSTRUMPACKSetHSSRelTol()`,
513*2ef1f0ffSBarry Smith           `MatSTRUMPACKSetHSSAbsTol()`, `MatSTRUMPACKSetHSSMaxRank()`, `MatSTRUMPACKSetHSSLeafSize()`, `MatSTRUMPACKSetHSSMinSepSize()`
514575704cbSPieter Ghysels M*/
515d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetFactor_aij_strumpack(Mat A, MatFactorType ftype, Mat *F)
516d71ae5a4SJacob Faibussowitsch {
5177d6ea485SPieter Ghysels   Mat                       B;
5187d6ea485SPieter Ghysels   PetscInt                  M = A->rmap->N, N = A->cmap->N;
51926cc229bSBarry Smith   PetscBool                 verb, flg;
52034f43fa5SPieter Ghysels   STRUMPACK_SparseSolver   *S;
52134f43fa5SPieter Ghysels   STRUMPACK_INTERFACE       iface;
5229371c9d4SSatish Balay   const STRUMPACK_PRECISION table[2][2][2] = {
5239371c9d4SSatish Balay     {{STRUMPACK_FLOATCOMPLEX_64, STRUMPACK_DOUBLECOMPLEX_64}, {STRUMPACK_FLOAT_64, STRUMPACK_DOUBLE_64}},
5249371c9d4SSatish Balay     {{STRUMPACK_FLOATCOMPLEX, STRUMPACK_DOUBLECOMPLEX},       {STRUMPACK_FLOAT, STRUMPACK_DOUBLE}      }
5259371c9d4SSatish Balay   };
5264ac6704cSBarry Smith   const STRUMPACK_PRECISION prec = table[(sizeof(PetscInt) == 8) ? 0 : 1][(PETSC_SCALAR == PETSC_COMPLEX) ? 0 : 1][(PETSC_REAL == PETSC_FLOAT) ? 0 : 1];
5277d6ea485SPieter Ghysels 
5287d6ea485SPieter Ghysels   PetscFunctionBegin;
5297d6ea485SPieter Ghysels   /* Create the factorization matrix */
5309566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
5319566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, M, N));
53235e8bcc9SJunchao Zhang   PetscCall(PetscStrallocpy("strumpack", &((PetscObject)B)->type_name));
53335e8bcc9SJunchao Zhang   PetscCall(MatSetUp(B));
53466e17bc3SBarry Smith   B->trivialsymbolic = PETSC_TRUE;
535575704cbSPieter Ghysels   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) {
5367d6ea485SPieter Ghysels     B->ops->lufactorsymbolic  = MatLUFactorSymbolic_STRUMPACK;
537575704cbSPieter Ghysels     B->ops->ilufactorsymbolic = MatLUFactorSymbolic_STRUMPACK;
538575704cbSPieter Ghysels   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Factor type not supported");
53935e8bcc9SJunchao Zhang   B->ops->getinfo     = MatGetInfo_External;
5407d6ea485SPieter Ghysels   B->ops->view        = MatView_STRUMPACK;
5417d6ea485SPieter Ghysels   B->ops->destroy     = MatDestroy_STRUMPACK;
5427d6ea485SPieter Ghysels   B->ops->getdiagonal = MatGetDiagonal_STRUMPACK;
5439566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_aij_strumpack));
5449566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetReordering_C", MatSTRUMPACKSetReordering_STRUMPACK));
5459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetColPerm_C", MatSTRUMPACKSetColPerm_STRUMPACK));
5469566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetHSSRelTol_C", MatSTRUMPACKSetHSSRelTol_STRUMPACK));
5479566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetHSSAbsTol_C", MatSTRUMPACKSetHSSAbsTol_STRUMPACK));
5489566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetHSSMaxRank_C", MatSTRUMPACKSetHSSMaxRank_STRUMPACK));
5499566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetHSSLeafSize_C", MatSTRUMPACKSetHSSLeafSize_STRUMPACK));
5509566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetHSSMinSepSize_C", MatSTRUMPACKSetHSSMinSepSize_STRUMPACK));
551575704cbSPieter Ghysels   B->factortype = ftype;
5524dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&S));
55335e8bcc9SJunchao Zhang   B->data = S;
5540d08a34dSPieter Ghysels 
55535e8bcc9SJunchao Zhang   PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATSEQAIJ, &flg)); /* A might be MATSEQAIJCUSPARSE */
5560d08a34dSPieter Ghysels   iface = flg ? STRUMPACK_MT : STRUMPACK_MPI_DIST;
5577d6ea485SPieter Ghysels 
55826cc229bSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)B), ((PetscObject)B)->prefix, "STRUMPACK Options", "Mat");
559575704cbSPieter Ghysels   verb = PetscLogPrintInfo ? PETSC_TRUE : PETSC_FALSE;
5609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_strumpack_verbose", "Print STRUMPACK information", "None", verb, &verb, NULL));
5617d6ea485SPieter Ghysels 
562e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_init", STRUMPACK_init(S, PetscObjectComm((PetscObject)A), prec, iface, 0, NULL, verb));
56355c022e5SPieter Ghysels 
56488382b05SPieter Ghysels   if (ftype == MAT_FACTOR_ILU) {
56588382b05SPieter Ghysels     /* When enabling HSS compression, the STRUMPACK solver becomes an incomplete              */
56688382b05SPieter Ghysels     /* (or approximate) LU factorization.                                                     */
567e77caa6dSBarry Smith     PetscStackCallExternalVoid("STRUMPACK_enable_HSS", STRUMPACK_enable_HSS(*S));
56888382b05SPieter Ghysels   }
569d0609cedSBarry Smith   PetscOptionsEnd();
57055c022e5SPieter Ghysels 
57126cc229bSBarry Smith   /* set solvertype */
57226cc229bSBarry Smith   PetscCall(PetscFree(B->solvertype));
57326cc229bSBarry Smith   PetscCall(PetscStrallocpy(MATSOLVERSTRUMPACK, &B->solvertype));
57426cc229bSBarry Smith 
5757d6ea485SPieter Ghysels   *F = B;
5763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5777d6ea485SPieter Ghysels }
5787d6ea485SPieter Ghysels 
579d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_STRUMPACK(void)
580d71ae5a4SJacob Faibussowitsch {
5817d6ea485SPieter Ghysels   PetscFunctionBegin;
5829566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERSTRUMPACK, MATMPIAIJ, MAT_FACTOR_LU, MatGetFactor_aij_strumpack));
5839566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERSTRUMPACK, MATSEQAIJ, MAT_FACTOR_LU, MatGetFactor_aij_strumpack));
5849566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERSTRUMPACK, MATMPIAIJ, MAT_FACTOR_ILU, MatGetFactor_aij_strumpack));
5859566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERSTRUMPACK, MATSEQAIJ, MAT_FACTOR_ILU, MatGetFactor_aij_strumpack));
5863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5877d6ea485SPieter Ghysels }
588