xref: /petsc/src/mat/impls/aij/mpi/strumpack/strumpack.c (revision c3339decea92175325d9368fa13196bcd0e0e58b)
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");
97d6ea485SPieter Ghysels   PetscFunctionReturn(0);
107d6ea485SPieter Ghysels }
117d6ea485SPieter Ghysels 
12d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_STRUMPACK(Mat A)
13d71ae5a4SJacob Faibussowitsch {
140d08a34dSPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)A->spptr;
157d6ea485SPieter Ghysels   PetscBool               flg;
167d6ea485SPieter Ghysels 
177d6ea485SPieter Ghysels   PetscFunctionBegin;
187d6ea485SPieter Ghysels   /* Deallocate STRUMPACK storage */
19e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_destroy", STRUMPACK_destroy(S));
209566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->spptr));
219566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &flg));
227d6ea485SPieter Ghysels   if (flg) {
239566063dSJacob Faibussowitsch     PetscCall(MatDestroy_SeqAIJ(A));
247d6ea485SPieter Ghysels   } else {
259566063dSJacob Faibussowitsch     PetscCall(MatDestroy_MPIAIJ(A));
267d6ea485SPieter Ghysels   }
27575704cbSPieter Ghysels 
28575704cbSPieter Ghysels   /* clear composed functions */
299566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorGetSolverType_C", NULL));
309566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetReordering_C", NULL));
319566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetColPerm_C", NULL));
329566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetHSSRelTol_C", NULL));
339566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetHSSAbsTol_C", NULL));
349566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetHSSMaxRank_C", NULL));
359566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetHSSLeafSize_C", NULL));
369566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSTRUMPACKSetHSSMinSepSize_C", NULL));
37575704cbSPieter Ghysels 
38575704cbSPieter Ghysels   PetscFunctionReturn(0);
39575704cbSPieter Ghysels }
40575704cbSPieter Ghysels 
41d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSTRUMPACKSetReordering_STRUMPACK(Mat F, MatSTRUMPACKReordering reordering)
42d71ae5a4SJacob Faibussowitsch {
430d08a34dSPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->spptr;
44575704cbSPieter Ghysels 
45575704cbSPieter Ghysels   PetscFunctionBegin;
46e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_reordering_method", STRUMPACK_set_reordering_method(*S, (STRUMPACK_REORDERING_STRATEGY)reordering));
4734f43fa5SPieter Ghysels   PetscFunctionReturn(0);
4834f43fa5SPieter Ghysels }
49e5a36eccSStefano Zampini 
50542aee0fSPieter Ghysels /*@
5134f43fa5SPieter Ghysels   MatSTRUMPACKSetReordering - Set STRUMPACK fill-reducing reordering
5234f43fa5SPieter Ghysels 
5334f43fa5SPieter Ghysels    Input Parameters:
5411a5261eSBarry Smith +  F - the factored matrix obtained by calling `MatGetFactor(`) from PETSc-STRUMPACK interface
5534f43fa5SPieter Ghysels -  reordering - the code to be used to find the fill-reducing reordering
5634f43fa5SPieter Ghysels       Possible values: NATURAL=0 METIS=1 PARMETIS=2 SCOTCH=3 PTSCOTCH=4 RCM=5
5734f43fa5SPieter Ghysels 
5811a5261eSBarry Smith   Options Database Key:
5934f43fa5SPieter Ghysels .   -mat_strumpack_reordering <METIS>  - Sparsity reducing matrix reordering (choose one of) NATURAL METIS PARMETIS SCOTCH PTSCOTCH RCM (None)
6034f43fa5SPieter Ghysels 
6134f43fa5SPieter Ghysels    Level: beginner
6234f43fa5SPieter Ghysels 
6334f43fa5SPieter Ghysels    References:
64606c0280SSatish Balay .  * - STRUMPACK manual
6534f43fa5SPieter Ghysels 
66db781477SPatrick Sanan .seealso: `MatGetFactor()`
67542aee0fSPieter Ghysels @*/
68d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSTRUMPACKSetReordering(Mat F, MatSTRUMPACKReordering reordering)
69d71ae5a4SJacob Faibussowitsch {
7034f43fa5SPieter Ghysels   PetscFunctionBegin;
7134f43fa5SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
7234f43fa5SPieter Ghysels   PetscValidLogicalCollectiveEnum(F, reordering, 2);
73cac4c232SBarry Smith   PetscTryMethod(F, "MatSTRUMPACKSetReordering_C", (Mat, MatSTRUMPACKReordering), (F, reordering));
74575704cbSPieter Ghysels   PetscFunctionReturn(0);
75575704cbSPieter Ghysels }
76575704cbSPieter Ghysels 
77d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSTRUMPACKSetColPerm_STRUMPACK(Mat F, PetscBool cperm)
78d71ae5a4SJacob Faibussowitsch {
7934f43fa5SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->spptr;
8034f43fa5SPieter Ghysels 
8134f43fa5SPieter Ghysels   PetscFunctionBegin;
82e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_set_mc64job", STRUMPACK_set_mc64job(*S, cperm ? 5 : 0));
8334f43fa5SPieter Ghysels   PetscFunctionReturn(0);
8434f43fa5SPieter Ghysels }
85e5a36eccSStefano Zampini 
86575704cbSPieter Ghysels /*@
8734f43fa5SPieter Ghysels   MatSTRUMPACKSetColPerm - Set whether STRUMPACK should try to permute the columns of the matrix in order to get a nonzero diagonal
88147403d9SBarry Smith 
89*c3339decSBarry Smith    Logically Collective
90575704cbSPieter Ghysels 
91575704cbSPieter Ghysels    Input Parameters:
9211a5261eSBarry Smith +  F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface
9311a5261eSBarry Smith -  cperm - `PETSC_TRUE` to permute (internally) the columns of the matrix
94575704cbSPieter Ghysels 
9511a5261eSBarry Smith   Options Database Key:
96147403d9SBarry Smith .   -mat_strumpack_colperm <cperm> - true to use the permutation
97575704cbSPieter Ghysels 
98575704cbSPieter Ghysels    Level: beginner
99575704cbSPieter Ghysels 
100575704cbSPieter Ghysels    References:
101606c0280SSatish Balay .  * - STRUMPACK manual
102575704cbSPieter Ghysels 
103db781477SPatrick Sanan .seealso: `MatGetFactor()`
104575704cbSPieter Ghysels @*/
105d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSTRUMPACKSetColPerm(Mat F, PetscBool cperm)
106d71ae5a4SJacob Faibussowitsch {
107575704cbSPieter Ghysels   PetscFunctionBegin;
108575704cbSPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
10934f43fa5SPieter Ghysels   PetscValidLogicalCollectiveBool(F, cperm, 2);
110cac4c232SBarry Smith   PetscTryMethod(F, "MatSTRUMPACKSetColPerm_C", (Mat, PetscBool), (F, cperm));
111575704cbSPieter Ghysels   PetscFunctionReturn(0);
112575704cbSPieter Ghysels }
113575704cbSPieter Ghysels 
114d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSTRUMPACKSetHSSRelTol_STRUMPACK(Mat F, PetscReal rtol)
115d71ae5a4SJacob Faibussowitsch {
11634f43fa5SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->spptr;
11734f43fa5SPieter Ghysels 
11834f43fa5SPieter Ghysels   PetscFunctionBegin;
119e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_set_HSS_rel_tol", STRUMPACK_set_HSS_rel_tol(*S, rtol));
12034f43fa5SPieter Ghysels   PetscFunctionReturn(0);
12134f43fa5SPieter Ghysels }
122e5a36eccSStefano Zampini 
12334f43fa5SPieter Ghysels /*@
12434f43fa5SPieter Ghysels   MatSTRUMPACKSetHSSRelTol - Set STRUMPACK relative tolerance for HSS compression
125147403d9SBarry Smith 
126*c3339decSBarry Smith   Logically Collective
12734f43fa5SPieter Ghysels 
12834f43fa5SPieter Ghysels    Input Parameters:
12911a5261eSBarry Smith +  F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface
13034f43fa5SPieter Ghysels -  rtol - relative compression tolerance
13134f43fa5SPieter Ghysels 
13211a5261eSBarry Smith   Options Database Key:
13334f43fa5SPieter Ghysels .   -mat_strumpack_hss_rel_tol <1e-2>         - Relative compression tolerance (None)
13434f43fa5SPieter Ghysels 
13534f43fa5SPieter Ghysels    Level: beginner
13634f43fa5SPieter Ghysels 
13734f43fa5SPieter Ghysels    References:
138606c0280SSatish Balay .  * - STRUMPACK manual
13934f43fa5SPieter Ghysels 
140db781477SPatrick Sanan .seealso: `MatGetFactor()`
14134f43fa5SPieter Ghysels @*/
142d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSTRUMPACKSetHSSRelTol(Mat F, PetscReal rtol)
143d71ae5a4SJacob Faibussowitsch {
14434f43fa5SPieter Ghysels   PetscFunctionBegin;
14534f43fa5SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
14634f43fa5SPieter Ghysels   PetscValidLogicalCollectiveReal(F, rtol, 2);
147cac4c232SBarry Smith   PetscTryMethod(F, "MatSTRUMPACKSetHSSRelTol_C", (Mat, PetscReal), (F, rtol));
14834f43fa5SPieter Ghysels   PetscFunctionReturn(0);
14934f43fa5SPieter Ghysels }
15034f43fa5SPieter Ghysels 
151d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSTRUMPACKSetHSSAbsTol_STRUMPACK(Mat F, PetscReal atol)
152d71ae5a4SJacob Faibussowitsch {
15334f43fa5SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->spptr;
15434f43fa5SPieter Ghysels 
15534f43fa5SPieter Ghysels   PetscFunctionBegin;
156e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_set_HSS_abs_tol", STRUMPACK_set_HSS_abs_tol(*S, atol));
15734f43fa5SPieter Ghysels   PetscFunctionReturn(0);
15834f43fa5SPieter Ghysels }
159e5a36eccSStefano Zampini 
16034f43fa5SPieter Ghysels /*@
16134f43fa5SPieter Ghysels   MatSTRUMPACKSetHSSAbsTol - Set STRUMPACK absolute tolerance for HSS compression
162147403d9SBarry Smith 
163*c3339decSBarry Smith    Logically Collective
16434f43fa5SPieter Ghysels 
16534f43fa5SPieter Ghysels    Input Parameters:
16611a5261eSBarry Smith +  F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface
16734f43fa5SPieter Ghysels -  atol - absolute compression tolerance
16834f43fa5SPieter Ghysels 
16911a5261eSBarry Smith   Options Database Key:
17034f43fa5SPieter Ghysels .   -mat_strumpack_hss_abs_tol <1e-8>         - Absolute compression tolerance (None)
17134f43fa5SPieter Ghysels 
17234f43fa5SPieter Ghysels    Level: beginner
17334f43fa5SPieter Ghysels 
17434f43fa5SPieter Ghysels    References:
175606c0280SSatish Balay .  * - STRUMPACK manual
17634f43fa5SPieter Ghysels 
177db781477SPatrick Sanan .seealso: `MatGetFactor()`
17834f43fa5SPieter Ghysels @*/
179d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSTRUMPACKSetHSSAbsTol(Mat F, PetscReal atol)
180d71ae5a4SJacob Faibussowitsch {
18134f43fa5SPieter Ghysels   PetscFunctionBegin;
18234f43fa5SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
18334f43fa5SPieter Ghysels   PetscValidLogicalCollectiveReal(F, atol, 2);
184cac4c232SBarry Smith   PetscTryMethod(F, "MatSTRUMPACKSetHSSAbsTol_C", (Mat, PetscReal), (F, atol));
18534f43fa5SPieter Ghysels   PetscFunctionReturn(0);
18634f43fa5SPieter Ghysels }
18734f43fa5SPieter Ghysels 
188d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSTRUMPACKSetHSSMaxRank_STRUMPACK(Mat F, PetscInt hssmaxrank)
189d71ae5a4SJacob Faibussowitsch {
19034f43fa5SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->spptr;
19134f43fa5SPieter Ghysels 
19234f43fa5SPieter Ghysels   PetscFunctionBegin;
193e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_set_HSS_max_rank", STRUMPACK_set_HSS_max_rank(*S, hssmaxrank));
19434f43fa5SPieter Ghysels   PetscFunctionReturn(0);
19534f43fa5SPieter Ghysels }
196e5a36eccSStefano Zampini 
19734f43fa5SPieter Ghysels /*@
19834f43fa5SPieter Ghysels   MatSTRUMPACKSetHSSMaxRank - Set STRUMPACK maximum HSS rank
199147403d9SBarry Smith 
200*c3339decSBarry Smith    Logically Collective
20134f43fa5SPieter Ghysels 
20234f43fa5SPieter Ghysels    Input Parameters:
20311a5261eSBarry Smith +  F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface
20434f43fa5SPieter Ghysels -  hssmaxrank - maximum rank used in low-rank approximation
20534f43fa5SPieter Ghysels 
20611a5261eSBarry Smith   Options Database Key:
20734f43fa5SPieter Ghysels .   -mat_strumpack_max_rank - Maximum rank in HSS compression, when using pctype ilu (None)
20834f43fa5SPieter Ghysels 
20934f43fa5SPieter Ghysels    Level: beginner
21034f43fa5SPieter Ghysels 
21134f43fa5SPieter Ghysels    References:
212606c0280SSatish Balay .  * - STRUMPACK manual
21334f43fa5SPieter Ghysels 
214db781477SPatrick Sanan .seealso: `MatGetFactor()`
21534f43fa5SPieter Ghysels @*/
216d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSTRUMPACKSetHSSMaxRank(Mat F, PetscInt hssmaxrank)
217d71ae5a4SJacob Faibussowitsch {
21834f43fa5SPieter Ghysels   PetscFunctionBegin;
21934f43fa5SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
22034f43fa5SPieter Ghysels   PetscValidLogicalCollectiveInt(F, hssmaxrank, 2);
221cac4c232SBarry Smith   PetscTryMethod(F, "MatSTRUMPACKSetHSSMaxRank_C", (Mat, PetscInt), (F, hssmaxrank));
22234f43fa5SPieter Ghysels   PetscFunctionReturn(0);
22334f43fa5SPieter Ghysels }
22434f43fa5SPieter Ghysels 
225d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSTRUMPACKSetHSSLeafSize_STRUMPACK(Mat F, PetscInt leaf_size)
226d71ae5a4SJacob Faibussowitsch {
227a36bf211SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->spptr;
228a36bf211SPieter Ghysels 
229a36bf211SPieter Ghysels   PetscFunctionBegin;
230e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_set_HSS_leaf_size", STRUMPACK_set_HSS_leaf_size(*S, leaf_size));
231a36bf211SPieter Ghysels   PetscFunctionReturn(0);
232a36bf211SPieter Ghysels }
233e5a36eccSStefano Zampini 
234a36bf211SPieter Ghysels /*@
235a36bf211SPieter Ghysels   MatSTRUMPACKSetHSSLeafSize - Set STRUMPACK HSS leaf size
236147403d9SBarry Smith 
237*c3339decSBarry Smith    Logically Collective
238a36bf211SPieter Ghysels 
239a36bf211SPieter Ghysels    Input Parameters:
24011a5261eSBarry Smith +  F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface
241a36bf211SPieter Ghysels -  leaf_size - Size of diagonal blocks in HSS approximation
242a36bf211SPieter Ghysels 
24311a5261eSBarry Smith   Options Database Key:
244a36bf211SPieter Ghysels .   -mat_strumpack_leaf_size - Size of diagonal blocks in HSS approximation, when using pctype ilu (None)
245a36bf211SPieter Ghysels 
246a36bf211SPieter Ghysels    Level: beginner
247a36bf211SPieter Ghysels 
248a36bf211SPieter Ghysels    References:
249606c0280SSatish Balay .  * - STRUMPACK manual
250a36bf211SPieter Ghysels 
251db781477SPatrick Sanan .seealso: `MatGetFactor()`
252a36bf211SPieter Ghysels @*/
253d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSTRUMPACKSetHSSLeafSize(Mat F, PetscInt leaf_size)
254d71ae5a4SJacob Faibussowitsch {
255a36bf211SPieter Ghysels   PetscFunctionBegin;
256a36bf211SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
257a36bf211SPieter Ghysels   PetscValidLogicalCollectiveInt(F, leaf_size, 2);
258cac4c232SBarry Smith   PetscTryMethod(F, "MatSTRUMPACKSetHSSLeafSize_C", (Mat, PetscInt), (F, leaf_size));
259a36bf211SPieter Ghysels   PetscFunctionReturn(0);
260a36bf211SPieter Ghysels }
261a36bf211SPieter Ghysels 
262d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSTRUMPACKSetHSSMinSepSize_STRUMPACK(Mat F, PetscInt hssminsize)
263d71ae5a4SJacob Faibussowitsch {
26434f43fa5SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->spptr;
26534f43fa5SPieter Ghysels 
26634f43fa5SPieter Ghysels   PetscFunctionBegin;
267e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_set_HSS_min_sep_size", STRUMPACK_set_HSS_min_sep_size(*S, hssminsize));
26834f43fa5SPieter Ghysels   PetscFunctionReturn(0);
26934f43fa5SPieter Ghysels }
270e5a36eccSStefano Zampini 
271291d0ed5SPieter Ghysels /*@
272291d0ed5SPieter Ghysels   MatSTRUMPACKSetHSSMinSepSize - Set STRUMPACK minimum separator size for low-rank approximation
273147403d9SBarry Smith 
274*c3339decSBarry Smith    Logically Collective
275291d0ed5SPieter Ghysels 
276291d0ed5SPieter Ghysels    Input Parameters:
27711a5261eSBarry Smith +  F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-STRUMPACK interface
278291d0ed5SPieter Ghysels -  hssminsize - minimum dense matrix size for low-rank approximation
279291d0ed5SPieter Ghysels 
28011a5261eSBarry Smith   Options Database Key:
281147403d9SBarry Smith .   -mat_strumpack_hss_min_sep_size <hssminsize> - set the minimum separator size
282291d0ed5SPieter Ghysels 
283291d0ed5SPieter Ghysels    Level: beginner
284291d0ed5SPieter Ghysels 
285291d0ed5SPieter Ghysels    References:
286606c0280SSatish Balay .  * - STRUMPACK manual
287291d0ed5SPieter Ghysels 
288db781477SPatrick Sanan .seealso: `MatGetFactor()`
289291d0ed5SPieter Ghysels @*/
290d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat F, PetscInt hssminsize)
291d71ae5a4SJacob Faibussowitsch {
292291d0ed5SPieter Ghysels   PetscFunctionBegin;
293291d0ed5SPieter Ghysels   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
294291d0ed5SPieter Ghysels   PetscValidLogicalCollectiveInt(F, hssminsize, 2);
295cac4c232SBarry Smith   PetscTryMethod(F, "MatSTRUMPACKSetHSSMinSepSize_C", (Mat, PetscInt), (F, hssminsize));
296291d0ed5SPieter Ghysels   PetscFunctionReturn(0);
297291d0ed5SPieter Ghysels }
298291d0ed5SPieter Ghysels 
299d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolve_STRUMPACK(Mat A, Vec b_mpi, Vec x)
300d71ae5a4SJacob Faibussowitsch {
3010d08a34dSPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)A->spptr;
3027d6ea485SPieter Ghysels   STRUMPACK_RETURN_CODE   sp_err;
3037d6ea485SPieter Ghysels   const PetscScalar      *bptr;
3047d6ea485SPieter Ghysels   PetscScalar            *xptr;
3057d6ea485SPieter Ghysels 
3067d6ea485SPieter Ghysels   PetscFunctionBegin;
3079566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x, &xptr));
3089566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(b_mpi, &bptr));
3090d08a34dSPieter Ghysels 
310e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_solve", sp_err = STRUMPACK_solve(*S, (PetscScalar *)bptr, xptr, 0));
3110d08a34dSPieter Ghysels   switch (sp_err) {
312d71ae5a4SJacob Faibussowitsch   case STRUMPACK_SUCCESS:
313d71ae5a4SJacob Faibussowitsch     break;
3149371c9d4SSatish Balay   case STRUMPACK_MATRIX_NOT_SET: {
3159371c9d4SSatish Balay     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix was not set");
3169371c9d4SSatish Balay     break;
3179371c9d4SSatish Balay   }
3189371c9d4SSatish Balay   case STRUMPACK_REORDERING_ERROR: {
3199371c9d4SSatish Balay     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix reordering failed");
3209371c9d4SSatish Balay     break;
3219371c9d4SSatish Balay   }
322d71ae5a4SJacob Faibussowitsch   default:
323d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: solve failed");
3247d6ea485SPieter Ghysels   }
3259566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x, &xptr));
3269566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(b_mpi, &bptr));
3277d6ea485SPieter Ghysels   PetscFunctionReturn(0);
3287d6ea485SPieter Ghysels }
3297d6ea485SPieter Ghysels 
330d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMatSolve_STRUMPACK(Mat A, Mat B_mpi, Mat X)
331d71ae5a4SJacob Faibussowitsch {
3327d6ea485SPieter Ghysels   PetscBool flg;
3337d6ea485SPieter Ghysels 
3347d6ea485SPieter Ghysels   PetscFunctionBegin;
3359566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)B_mpi, &flg, MATSEQDENSE, MATMPIDENSE, NULL));
3365f80ce2aSJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Matrix B must be MATDENSE matrix");
3379566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)X, &flg, MATSEQDENSE, MATMPIDENSE, NULL));
3385f80ce2aSJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Matrix X must be MATDENSE matrix");
3397d6ea485SPieter Ghysels   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "MatMatSolve_STRUMPACK() is not implemented yet");
3407d6ea485SPieter Ghysels   PetscFunctionReturn(0);
3417d6ea485SPieter Ghysels }
3427d6ea485SPieter Ghysels 
343d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_Info_STRUMPACK(Mat A, PetscViewer viewer)
344d71ae5a4SJacob Faibussowitsch {
345ad0c5e61SPieter Ghysels   PetscFunctionBegin;
346ad0c5e61SPieter Ghysels   /* check if matrix is strumpack type */
347ad0c5e61SPieter Ghysels   if (A->ops->solve != MatSolve_STRUMPACK) PetscFunctionReturn(0);
3489566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "STRUMPACK sparse solver!\n"));
349ad0c5e61SPieter Ghysels   PetscFunctionReturn(0);
350ad0c5e61SPieter Ghysels }
351ad0c5e61SPieter Ghysels 
352d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_STRUMPACK(Mat A, PetscViewer viewer)
353d71ae5a4SJacob Faibussowitsch {
354ad0c5e61SPieter Ghysels   PetscBool         iascii;
355ad0c5e61SPieter Ghysels   PetscViewerFormat format;
356ad0c5e61SPieter Ghysels 
357ad0c5e61SPieter Ghysels   PetscFunctionBegin;
3589566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
359ad0c5e61SPieter Ghysels   if (iascii) {
3609566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
3619566063dSJacob Faibussowitsch     if (format == PETSC_VIEWER_ASCII_INFO) PetscCall(MatView_Info_STRUMPACK(A, viewer));
362ad0c5e61SPieter Ghysels   }
363ad0c5e61SPieter Ghysels   PetscFunctionReturn(0);
364ad0c5e61SPieter Ghysels }
3657d6ea485SPieter Ghysels 
366d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLUFactorNumeric_STRUMPACK(Mat F, Mat A, const MatFactorInfo *info)
367d71ae5a4SJacob Faibussowitsch {
3680d08a34dSPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver *)F->spptr;
3697d6ea485SPieter Ghysels   STRUMPACK_RETURN_CODE   sp_err;
3707d6ea485SPieter Ghysels   Mat_SeqAIJ             *A_d, *A_o;
3717d6ea485SPieter Ghysels   Mat_MPIAIJ             *mat;
3720d08a34dSPieter Ghysels   PetscInt                M = A->rmap->N, m = A->rmap->n;
3737d6ea485SPieter Ghysels   PetscBool               flg;
3747d6ea485SPieter Ghysels 
3757d6ea485SPieter Ghysels   PetscFunctionBegin;
3769566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIAIJ, &flg));
3777d6ea485SPieter Ghysels   if (flg) { /* A is MATMPIAIJ */
3787d6ea485SPieter Ghysels     mat = (Mat_MPIAIJ *)A->data;
3797d6ea485SPieter Ghysels     A_d = (Mat_SeqAIJ *)(mat->A)->data;
3807d6ea485SPieter Ghysels     A_o = (Mat_SeqAIJ *)(mat->B)->data;
381e77caa6dSBarry Smith     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));
3827d6ea485SPieter Ghysels   } else { /* A is MATSEQAIJ */
3837d6ea485SPieter Ghysels     A_d = (Mat_SeqAIJ *)A->data;
384e77caa6dSBarry Smith     PetscStackCallExternalVoid("STRUMPACK_set_csr_matrix", STRUMPACK_set_csr_matrix(*S, &M, A_d->i, A_d->j, A_d->a, 0));
3857d6ea485SPieter Ghysels   }
3867d6ea485SPieter Ghysels 
3877d6ea485SPieter Ghysels   /* Reorder and Factor the matrix. */
3887d6ea485SPieter Ghysels   /* TODO figure out how to avoid reorder if the matrix values changed, but the pattern remains the same. */
389e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_reorder", sp_err = STRUMPACK_reorder(*S));
390e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_factor", sp_err = STRUMPACK_factor(*S));
3910d08a34dSPieter Ghysels   switch (sp_err) {
392d71ae5a4SJacob Faibussowitsch   case STRUMPACK_SUCCESS:
393d71ae5a4SJacob Faibussowitsch     break;
3949371c9d4SSatish Balay   case STRUMPACK_MATRIX_NOT_SET: {
3959371c9d4SSatish Balay     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix was not set");
3969371c9d4SSatish Balay     break;
3979371c9d4SSatish Balay   }
3989371c9d4SSatish Balay   case STRUMPACK_REORDERING_ERROR: {
3999371c9d4SSatish Balay     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: matrix reordering failed");
4009371c9d4SSatish Balay     break;
4019371c9d4SSatish Balay   }
402d71ae5a4SJacob Faibussowitsch   default:
403d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "STRUMPACK error: factorization failed");
4047d6ea485SPieter Ghysels   }
405cb250fa3SPieter Ghysels   F->assembled    = PETSC_TRUE;
406cb250fa3SPieter Ghysels   F->preallocated = PETSC_TRUE;
4077d6ea485SPieter Ghysels   PetscFunctionReturn(0);
4087d6ea485SPieter Ghysels }
4097d6ea485SPieter Ghysels 
410d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLUFactorSymbolic_STRUMPACK(Mat F, Mat A, IS r, IS c, const MatFactorInfo *info)
411d71ae5a4SJacob Faibussowitsch {
41226cc229bSBarry Smith   STRUMPACK_SparseSolver       *S = (STRUMPACK_SparseSolver *)F->spptr;
41326cc229bSBarry Smith   PetscBool                     flg, set;
41426cc229bSBarry Smith   PetscReal                     ctol;
41526cc229bSBarry Smith   PetscInt                      hssminsize, max_rank, leaf_size;
41626cc229bSBarry Smith   STRUMPACK_REORDERING_STRATEGY ndcurrent, ndvalue;
41726cc229bSBarry Smith   STRUMPACK_KRYLOV_SOLVER       itcurrent, itsolver;
41826cc229bSBarry Smith   const char *const             STRUMPACKNDTypes[] = {"NATURAL", "METIS", "PARMETIS", "SCOTCH", "PTSCOTCH", "RCM", "STRUMPACKNDTypes", "", 0};
41926cc229bSBarry Smith   const char *const             SolverTypes[]      = {"AUTO", "NONE", "REFINE", "PREC_GMRES", "GMRES", "PREC_BICGSTAB", "BICGSTAB", "SolverTypes", "", 0};
42026cc229bSBarry Smith 
4217d6ea485SPieter Ghysels   PetscFunctionBegin;
42226cc229bSBarry Smith   /* Set options to F */
42326cc229bSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)F), ((PetscObject)F)->prefix, "STRUMPACK Options", "Mat");
424e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_HSS_rel_tol", ctol = (PetscReal)STRUMPACK_HSS_rel_tol(*S));
42526cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_strumpack_hss_rel_tol", "Relative compression tolerance", "None", ctol, &ctol, &set));
426e77caa6dSBarry Smith   if (set) PetscStackCallExternalVoid("STRUMPACK_set_HSS_rel_tol", STRUMPACK_set_HSS_rel_tol(*S, (double)ctol));
42726cc229bSBarry Smith 
428e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_HSS_abs_tol", ctol = (PetscReal)STRUMPACK_HSS_abs_tol(*S));
42926cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_strumpack_hss_abs_tol", "Absolute compression tolerance", "None", ctol, &ctol, &set));
430e77caa6dSBarry Smith   if (set) PetscStackCallExternalVoid("STRUMPACK_set_HSS_abs_tol", STRUMPACK_set_HSS_abs_tol(*S, (double)ctol));
43126cc229bSBarry Smith 
432e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_mc64job", flg = (STRUMPACK_mc64job(*S) == 0) ? PETSC_FALSE : PETSC_TRUE);
43326cc229bSBarry Smith   PetscCall(PetscOptionsBool("-mat_strumpack_colperm", "Find a col perm to get nonzero diagonal", "None", flg, &flg, &set));
434e77caa6dSBarry Smith   if (set) PetscStackCallExternalVoid("STRUMPACK_set_mc64job", STRUMPACK_set_mc64job(*S, flg ? 5 : 0));
43526cc229bSBarry Smith 
436e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_HSS_min_sep_size", hssminsize = (PetscInt)STRUMPACK_HSS_min_sep_size(*S));
43726cc229bSBarry Smith   PetscCall(PetscOptionsInt("-mat_strumpack_hss_min_sep_size", "Minimum size of separator for HSS compression", "None", hssminsize, &hssminsize, &set));
438e77caa6dSBarry Smith   if (set) PetscStackCallExternalVoid("STRUMPACK_set_HSS_min_sep_size", STRUMPACK_set_HSS_min_sep_size(*S, (int)hssminsize));
43926cc229bSBarry Smith 
440e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_HSS_max_rank", max_rank = (PetscInt)STRUMPACK_HSS_max_rank(*S));
44126cc229bSBarry Smith   PetscCall(PetscOptionsInt("-mat_strumpack_max_rank", "Maximum rank in HSS approximation", "None", max_rank, &max_rank, &set));
442e77caa6dSBarry Smith   if (set) PetscStackCallExternalVoid("STRUMPACK_set_HSS_max_rank", STRUMPACK_set_HSS_max_rank(*S, (int)max_rank));
44326cc229bSBarry Smith 
444e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_HSS_leaf_size", leaf_size = (PetscInt)STRUMPACK_HSS_leaf_size(*S));
44526cc229bSBarry Smith   PetscCall(PetscOptionsInt("-mat_strumpack_leaf_size", "Size of diagonal blocks in HSS approximation", "None", leaf_size, &leaf_size, &set));
446e77caa6dSBarry Smith   if (set) PetscStackCallExternalVoid("STRUMPACK_set_HSS_leaf_size", STRUMPACK_set_HSS_leaf_size(*S, (int)leaf_size));
44726cc229bSBarry Smith 
448e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_reordering_method", ndcurrent = STRUMPACK_reordering_method(*S));
44926cc229bSBarry Smith   PetscCall(PetscOptionsEnum("-mat_strumpack_reordering", "Sparsity reducing matrix reordering", "None", STRUMPACKNDTypes, (PetscEnum)ndcurrent, (PetscEnum *)&ndvalue, &set));
450e77caa6dSBarry Smith   if (set) PetscStackCallExternalVoid("STRUMPACK_set_reordering_method", STRUMPACK_set_reordering_method(*S, ndvalue));
45126cc229bSBarry Smith 
45226cc229bSBarry Smith   /* Disable the outer iterative solver from STRUMPACK.                                       */
45326cc229bSBarry Smith   /* When STRUMPACK is used as a direct solver, it will by default do iterative refinement.   */
45426cc229bSBarry Smith   /* When STRUMPACK is used as an approximate factorization preconditioner (by enabling       */
45526cc229bSBarry Smith   /* low-rank compression), it will use it's own preconditioned GMRES. Here we can disable    */
45626cc229bSBarry Smith   /* the outer iterative solver, as PETSc uses STRUMPACK from within a KSP.                   */
457e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_set_Krylov_solver", STRUMPACK_set_Krylov_solver(*S, STRUMPACK_DIRECT));
45826cc229bSBarry Smith 
459e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_Krylov_solver", itcurrent = STRUMPACK_Krylov_solver(*S));
46026cc229bSBarry Smith   PetscCall(PetscOptionsEnum("-mat_strumpack_iterative_solver", "Select iterative solver from STRUMPACK", "None", SolverTypes, (PetscEnum)itcurrent, (PetscEnum *)&itsolver, &set));
461e77caa6dSBarry Smith   if (set) PetscStackCallExternalVoid("STRUMPACK_set_Krylov_solver", STRUMPACK_set_Krylov_solver(*S, itsolver));
46226cc229bSBarry Smith   PetscOptionsEnd();
46326cc229bSBarry Smith 
4647d6ea485SPieter Ghysels   F->ops->lufactornumeric = MatLUFactorNumeric_STRUMPACK;
4657d6ea485SPieter Ghysels   F->ops->solve           = MatSolve_STRUMPACK;
4667d6ea485SPieter Ghysels   F->ops->matsolve        = MatMatSolve_STRUMPACK;
4677d6ea485SPieter Ghysels   PetscFunctionReturn(0);
4687d6ea485SPieter Ghysels }
4697d6ea485SPieter Ghysels 
470d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorGetSolverType_aij_strumpack(Mat A, MatSolverType *type)
471d71ae5a4SJacob Faibussowitsch {
4727d6ea485SPieter Ghysels   PetscFunctionBegin;
4737d6ea485SPieter Ghysels   *type = MATSOLVERSTRUMPACK;
4747d6ea485SPieter Ghysels   PetscFunctionReturn(0);
4757d6ea485SPieter Ghysels }
4767d6ea485SPieter Ghysels 
477575704cbSPieter Ghysels /*MC
478575704cbSPieter Ghysels   MATSOLVERSSTRUMPACK = "strumpack" - A solver package providing a direct sparse solver (PCLU)
479575704cbSPieter Ghysels   and a preconditioner (PCILU) using low-rank compression via the external package STRUMPACK.
480575704cbSPieter Ghysels 
481575704cbSPieter Ghysels   Consult the STRUMPACK-sparse manual for more info.
482575704cbSPieter Ghysels 
483575704cbSPieter Ghysels   Use
484575704cbSPieter Ghysels      ./configure --download-strumpack
485575704cbSPieter Ghysels   to have PETSc installed with STRUMPACK
486575704cbSPieter Ghysels 
487575704cbSPieter Ghysels   Use
4883ca39a21SBarry Smith     -pc_type lu -pc_factor_mat_solver_type strumpack
489575704cbSPieter Ghysels   to use this as an exact (direct) solver, use
4903ca39a21SBarry Smith     -pc_type ilu -pc_factor_mat_solver_type strumpack
491575704cbSPieter Ghysels   to enable low-rank compression (i.e, use as a preconditioner).
492575704cbSPieter Ghysels 
493575704cbSPieter Ghysels   Works with AIJ matrices
494575704cbSPieter Ghysels 
495575704cbSPieter Ghysels   Options Database Keys:
49667b8a455SSatish Balay + -mat_strumpack_verbose                    - verbose info
49734f43fa5SPieter Ghysels . -mat_strumpack_hss_rel_tol <1e-2>         - Relative compression tolerance (None)
49834f43fa5SPieter Ghysels . -mat_strumpack_hss_abs_tol <1e-8>         - Absolute compression tolerance (None)
499575704cbSPieter Ghysels . -mat_strumpack_colperm <TRUE>             - Permute matrix to make diagonal nonzeros (None)
50034f43fa5SPieter Ghysels . -mat_strumpack_hss_min_sep_size <256>     - Minimum size of separator for HSS compression (None)
50134f43fa5SPieter Ghysels . -mat_strumpack_max_rank                   - Maximum rank in HSS compression, when using pctype ilu (None)
502a36bf211SPieter Ghysels . -mat_strumpack_leaf_size                  - Size of diagonal blocks in HSS approximation, when using pctype ilu (None)
50334f43fa5SPieter Ghysels . -mat_strumpack_reordering <METIS>         - Sparsity reducing matrix reordering (choose one of) NATURAL METIS PARMETIS SCOTCH PTSCOTCH RCM (None)
504a2b725a8SWilliam Gropp - -mat_strumpack_iterative_solver <DIRECT>  - Select iterative solver from STRUMPACK (choose one of) AUTO DIRECT REFINE PREC_GMRES GMRES PREC_BICGSTAB BICGSTAB (None)
505575704cbSPieter Ghysels 
506575704cbSPieter Ghysels  Level: beginner
507575704cbSPieter Ghysels 
50811a5261eSBarry Smith .seealso: `PCLU`, `PCILU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`, `PCFactorSetMatSolverType()`, `MatSolverType`, `MatGetFactor()`
509575704cbSPieter Ghysels M*/
510d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetFactor_aij_strumpack(Mat A, MatFactorType ftype, Mat *F)
511d71ae5a4SJacob Faibussowitsch {
5127d6ea485SPieter Ghysels   Mat                       B;
5137d6ea485SPieter Ghysels   PetscInt                  M = A->rmap->N, N = A->cmap->N;
51426cc229bSBarry Smith   PetscBool                 verb, flg;
51534f43fa5SPieter Ghysels   STRUMPACK_SparseSolver   *S;
51634f43fa5SPieter Ghysels   STRUMPACK_INTERFACE       iface;
5179371c9d4SSatish Balay   const STRUMPACK_PRECISION table[2][2][2] = {
5189371c9d4SSatish Balay     {{STRUMPACK_FLOATCOMPLEX_64, STRUMPACK_DOUBLECOMPLEX_64}, {STRUMPACK_FLOAT_64, STRUMPACK_DOUBLE_64}},
5199371c9d4SSatish Balay     {{STRUMPACK_FLOATCOMPLEX, STRUMPACK_DOUBLECOMPLEX},       {STRUMPACK_FLOAT, STRUMPACK_DOUBLE}      }
5209371c9d4SSatish Balay   };
5214ac6704cSBarry Smith   const STRUMPACK_PRECISION prec = table[(sizeof(PetscInt) == 8) ? 0 : 1][(PETSC_SCALAR == PETSC_COMPLEX) ? 0 : 1][(PETSC_REAL == PETSC_FLOAT) ? 0 : 1];
5227d6ea485SPieter Ghysels 
5237d6ea485SPieter Ghysels   PetscFunctionBegin;
5247d6ea485SPieter Ghysels   /* Create the factorization matrix */
5259566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
5269566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, M, N));
5279566063dSJacob Faibussowitsch   PetscCall(MatSetType(B, ((PetscObject)A)->type_name));
5289566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(B, 0, NULL));
5299566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(B, 0, NULL, 0, NULL));
53066e17bc3SBarry Smith   B->trivialsymbolic = PETSC_TRUE;
531575704cbSPieter Ghysels   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) {
5327d6ea485SPieter Ghysels     B->ops->lufactorsymbolic  = MatLUFactorSymbolic_STRUMPACK;
533575704cbSPieter Ghysels     B->ops->ilufactorsymbolic = MatLUFactorSymbolic_STRUMPACK;
534575704cbSPieter Ghysels   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Factor type not supported");
5357d6ea485SPieter Ghysels   B->ops->view        = MatView_STRUMPACK;
5367d6ea485SPieter Ghysels   B->ops->destroy     = MatDestroy_STRUMPACK;
5377d6ea485SPieter Ghysels   B->ops->getdiagonal = MatGetDiagonal_STRUMPACK;
5389566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_aij_strumpack));
5399566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetReordering_C", MatSTRUMPACKSetReordering_STRUMPACK));
5409566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetColPerm_C", MatSTRUMPACKSetColPerm_STRUMPACK));
5419566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetHSSRelTol_C", MatSTRUMPACKSetHSSRelTol_STRUMPACK));
5429566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetHSSAbsTol_C", MatSTRUMPACKSetHSSAbsTol_STRUMPACK));
5439566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetHSSMaxRank_C", MatSTRUMPACKSetHSSMaxRank_STRUMPACK));
5449566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetHSSLeafSize_C", MatSTRUMPACKSetHSSLeafSize_STRUMPACK));
5459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSTRUMPACKSetHSSMinSepSize_C", MatSTRUMPACKSetHSSMinSepSize_STRUMPACK));
546575704cbSPieter Ghysels   B->factortype = ftype;
5474dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&S));
5480d08a34dSPieter Ghysels   B->spptr = S;
5490d08a34dSPieter Ghysels 
5509566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &flg));
5510d08a34dSPieter Ghysels   iface = flg ? STRUMPACK_MT : STRUMPACK_MPI_DIST;
5527d6ea485SPieter Ghysels 
55326cc229bSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)B), ((PetscObject)B)->prefix, "STRUMPACK Options", "Mat");
554575704cbSPieter Ghysels   verb = PetscLogPrintInfo ? PETSC_TRUE : PETSC_FALSE;
5559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_strumpack_verbose", "Print STRUMPACK information", "None", verb, &verb, NULL));
5567d6ea485SPieter Ghysels 
557e77caa6dSBarry Smith   PetscStackCallExternalVoid("STRUMPACK_init", STRUMPACK_init(S, PetscObjectComm((PetscObject)A), prec, iface, 0, NULL, verb));
55855c022e5SPieter Ghysels 
55988382b05SPieter Ghysels   if (ftype == MAT_FACTOR_ILU) {
56088382b05SPieter Ghysels     /* When enabling HSS compression, the STRUMPACK solver becomes an incomplete              */
56188382b05SPieter Ghysels     /* (or approximate) LU factorization.                                                     */
562e77caa6dSBarry Smith     PetscStackCallExternalVoid("STRUMPACK_enable_HSS", STRUMPACK_enable_HSS(*S));
56388382b05SPieter Ghysels   }
564d0609cedSBarry Smith   PetscOptionsEnd();
56555c022e5SPieter Ghysels 
56626cc229bSBarry Smith   /* set solvertype */
56726cc229bSBarry Smith   PetscCall(PetscFree(B->solvertype));
56826cc229bSBarry Smith   PetscCall(PetscStrallocpy(MATSOLVERSTRUMPACK, &B->solvertype));
56926cc229bSBarry Smith 
5707d6ea485SPieter Ghysels   *F = B;
5717d6ea485SPieter Ghysels   PetscFunctionReturn(0);
5727d6ea485SPieter Ghysels }
5737d6ea485SPieter Ghysels 
574d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_STRUMPACK(void)
575d71ae5a4SJacob Faibussowitsch {
5767d6ea485SPieter Ghysels   PetscFunctionBegin;
5779566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERSTRUMPACK, MATMPIAIJ, MAT_FACTOR_LU, MatGetFactor_aij_strumpack));
5789566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERSTRUMPACK, MATSEQAIJ, MAT_FACTOR_LU, MatGetFactor_aij_strumpack));
5799566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERSTRUMPACK, MATMPIAIJ, MAT_FACTOR_ILU, MatGetFactor_aij_strumpack));
5809566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERSTRUMPACK, MATSEQAIJ, MAT_FACTOR_ILU, MatGetFactor_aij_strumpack));
5817d6ea485SPieter Ghysels   PetscFunctionReturn(0);
5827d6ea485SPieter Ghysels }
583