xref: /petsc/src/mat/impls/aij/mpi/strumpack/strumpack.c (revision 147403d98689287189d69e47992b7c8152b2c9da)
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 
5ad0c5e61SPieter Ghysels static PetscErrorCode MatGetDiagonal_STRUMPACK(Mat A,Vec v)
67d6ea485SPieter Ghysels {
77d6ea485SPieter Ghysels   PetscFunctionBegin;
87d6ea485SPieter Ghysels   SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type: STRUMPACK factor");
97d6ea485SPieter Ghysels   PetscFunctionReturn(0);
107d6ea485SPieter Ghysels }
117d6ea485SPieter Ghysels 
12ad0c5e61SPieter Ghysels static PetscErrorCode MatDestroy_STRUMPACK(Mat A)
137d6ea485SPieter Ghysels {
140d08a34dSPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)A->spptr;
157d6ea485SPieter Ghysels   PetscErrorCode         ierr;
167d6ea485SPieter Ghysels   PetscBool              flg;
177d6ea485SPieter Ghysels 
187d6ea485SPieter Ghysels   PetscFunctionBegin;
197d6ea485SPieter Ghysels   /* Deallocate STRUMPACK storage */
200d08a34dSPieter Ghysels   PetscStackCall("STRUMPACK_destroy",STRUMPACK_destroy(S));
217d6ea485SPieter Ghysels   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
227d6ea485SPieter Ghysels   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr);
237d6ea485SPieter Ghysels   if (flg) {
247d6ea485SPieter Ghysels     ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
257d6ea485SPieter Ghysels   } else {
267d6ea485SPieter Ghysels     ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr);
277d6ea485SPieter Ghysels   }
28575704cbSPieter Ghysels 
29575704cbSPieter Ghysels   /* clear composed functions */
303ca39a21SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverType_C",NULL);CHKERRQ(ierr);
3134f43fa5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetReordering_C",NULL);CHKERRQ(ierr);
3234f43fa5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetColPerm_C",NULL);CHKERRQ(ierr);
3334f43fa5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSRelTol_C",NULL);CHKERRQ(ierr);
3434f43fa5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSAbsTol_C",NULL);CHKERRQ(ierr);
3534f43fa5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSMaxRank_C",NULL);CHKERRQ(ierr);
36a36bf211SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSLeafSize_C",NULL);CHKERRQ(ierr);
37291d0ed5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSMinSepSize_C",NULL);CHKERRQ(ierr);
38575704cbSPieter Ghysels 
39575704cbSPieter Ghysels   PetscFunctionReturn(0);
40575704cbSPieter Ghysels }
41575704cbSPieter Ghysels 
4234f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetReordering_STRUMPACK(Mat F,MatSTRUMPACKReordering reordering)
43575704cbSPieter Ghysels {
440d08a34dSPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
45575704cbSPieter Ghysels 
46575704cbSPieter Ghysels   PetscFunctionBegin;
4734f43fa5SPieter Ghysels   PetscStackCall("STRUMPACK_reordering_method",STRUMPACK_set_reordering_method(*S,(STRUMPACK_REORDERING_STRATEGY)reordering));
4834f43fa5SPieter Ghysels   PetscFunctionReturn(0);
4934f43fa5SPieter Ghysels }
50e5a36eccSStefano Zampini 
51542aee0fSPieter Ghysels /*@
5234f43fa5SPieter Ghysels   MatSTRUMPACKSetReordering - Set STRUMPACK fill-reducing reordering
5334f43fa5SPieter Ghysels 
5434f43fa5SPieter Ghysels    Input Parameters:
5534f43fa5SPieter Ghysels +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
5634f43fa5SPieter Ghysels -  reordering - the code to be used to find the fill-reducing reordering
5734f43fa5SPieter Ghysels       Possible values: NATURAL=0 METIS=1 PARMETIS=2 SCOTCH=3 PTSCOTCH=4 RCM=5
5834f43fa5SPieter Ghysels 
5934f43fa5SPieter Ghysels   Options Database:
6034f43fa5SPieter Ghysels .   -mat_strumpack_reordering <METIS>  - Sparsity reducing matrix reordering (choose one of) NATURAL METIS PARMETIS SCOTCH PTSCOTCH RCM (None)
6134f43fa5SPieter Ghysels 
6234f43fa5SPieter Ghysels    Level: beginner
6334f43fa5SPieter Ghysels 
6434f43fa5SPieter Ghysels    References:
65606c0280SSatish Balay .  * - STRUMPACK manual
6634f43fa5SPieter Ghysels 
6734f43fa5SPieter Ghysels .seealso: MatGetFactor()
68542aee0fSPieter Ghysels @*/
6934f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetReordering(Mat F,MatSTRUMPACKReordering reordering)
7034f43fa5SPieter Ghysels {
7134f43fa5SPieter Ghysels   PetscErrorCode ierr;
7234f43fa5SPieter Ghysels 
7334f43fa5SPieter Ghysels   PetscFunctionBegin;
7434f43fa5SPieter Ghysels   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
7534f43fa5SPieter Ghysels   PetscValidLogicalCollectiveEnum(F,reordering,2);
7634f43fa5SPieter Ghysels   ierr = PetscTryMethod(F,"MatSTRUMPACKSetReordering_C",(Mat,MatSTRUMPACKReordering),(F,reordering));CHKERRQ(ierr);
77575704cbSPieter Ghysels   PetscFunctionReturn(0);
78575704cbSPieter Ghysels }
79575704cbSPieter Ghysels 
8034f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetColPerm_STRUMPACK(Mat F,PetscBool cperm)
8134f43fa5SPieter Ghysels {
8234f43fa5SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
8334f43fa5SPieter Ghysels 
8434f43fa5SPieter Ghysels   PetscFunctionBegin;
8534f43fa5SPieter Ghysels   PetscStackCall("STRUMPACK_set_mc64job",STRUMPACK_set_mc64job(*S,cperm ? 5 : 0));
8634f43fa5SPieter Ghysels   PetscFunctionReturn(0);
8734f43fa5SPieter Ghysels }
88e5a36eccSStefano Zampini 
89575704cbSPieter Ghysels /*@
9034f43fa5SPieter Ghysels   MatSTRUMPACKSetColPerm - Set whether STRUMPACK should try to permute the columns of the matrix in order to get a nonzero diagonal
91*147403d9SBarry Smith 
92575704cbSPieter Ghysels    Logically Collective on Mat
93575704cbSPieter Ghysels 
94575704cbSPieter Ghysels    Input Parameters:
95575704cbSPieter Ghysels +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
96*147403d9SBarry Smith -  cperm - PETSC_TRUE to permute (internally) the columns of the matrix
97575704cbSPieter Ghysels 
98575704cbSPieter Ghysels   Options Database:
99*147403d9SBarry Smith .   -mat_strumpack_colperm <cperm> - true to use the permutation
100575704cbSPieter Ghysels 
101575704cbSPieter Ghysels    Level: beginner
102575704cbSPieter Ghysels 
103575704cbSPieter Ghysels    References:
104606c0280SSatish Balay .  * - STRUMPACK manual
105575704cbSPieter Ghysels 
106575704cbSPieter Ghysels .seealso: MatGetFactor()
107575704cbSPieter Ghysels @*/
10834f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetColPerm(Mat F,PetscBool cperm)
109575704cbSPieter Ghysels {
110575704cbSPieter Ghysels   PetscErrorCode ierr;
111575704cbSPieter Ghysels 
112575704cbSPieter Ghysels   PetscFunctionBegin;
113575704cbSPieter Ghysels   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
11434f43fa5SPieter Ghysels   PetscValidLogicalCollectiveBool(F,cperm,2);
11534f43fa5SPieter Ghysels   ierr = PetscTryMethod(F,"MatSTRUMPACKSetColPerm_C",(Mat,PetscBool),(F,cperm));CHKERRQ(ierr);
116575704cbSPieter Ghysels   PetscFunctionReturn(0);
117575704cbSPieter Ghysels }
118575704cbSPieter Ghysels 
11934f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSRelTol_STRUMPACK(Mat F,PetscReal rtol)
12034f43fa5SPieter Ghysels {
12134f43fa5SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
12234f43fa5SPieter Ghysels 
12334f43fa5SPieter Ghysels   PetscFunctionBegin;
124a36bf211SPieter Ghysels   PetscStackCall("STRUMPACK_set_HSS_rel_tol", STRUMPACK_set_HSS_rel_tol(*S,rtol));
12534f43fa5SPieter Ghysels   PetscFunctionReturn(0);
12634f43fa5SPieter Ghysels }
127e5a36eccSStefano Zampini 
12834f43fa5SPieter Ghysels /*@
12934f43fa5SPieter Ghysels   MatSTRUMPACKSetHSSRelTol - Set STRUMPACK relative tolerance for HSS compression
130*147403d9SBarry Smith 
13134f43fa5SPieter Ghysels   Logically Collective on Mat
13234f43fa5SPieter Ghysels 
13334f43fa5SPieter Ghysels    Input Parameters:
13434f43fa5SPieter Ghysels +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
13534f43fa5SPieter Ghysels -  rtol - relative compression tolerance
13634f43fa5SPieter Ghysels 
13734f43fa5SPieter Ghysels   Options Database:
13834f43fa5SPieter Ghysels .   -mat_strumpack_hss_rel_tol <1e-2>         - Relative compression tolerance (None)
13934f43fa5SPieter Ghysels 
14034f43fa5SPieter Ghysels    Level: beginner
14134f43fa5SPieter Ghysels 
14234f43fa5SPieter Ghysels    References:
143606c0280SSatish Balay .  * - STRUMPACK manual
14434f43fa5SPieter Ghysels 
14534f43fa5SPieter Ghysels .seealso: MatGetFactor()
14634f43fa5SPieter Ghysels @*/
14734f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSRelTol(Mat F,PetscReal rtol)
14834f43fa5SPieter Ghysels {
14934f43fa5SPieter Ghysels   PetscErrorCode ierr;
15034f43fa5SPieter Ghysels 
15134f43fa5SPieter Ghysels   PetscFunctionBegin;
15234f43fa5SPieter Ghysels   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
15334f43fa5SPieter Ghysels   PetscValidLogicalCollectiveReal(F,rtol,2);
15434f43fa5SPieter Ghysels   ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSRelTol_C",(Mat,PetscReal),(F,rtol));CHKERRQ(ierr);
15534f43fa5SPieter Ghysels   PetscFunctionReturn(0);
15634f43fa5SPieter Ghysels }
15734f43fa5SPieter Ghysels 
15834f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSAbsTol_STRUMPACK(Mat F,PetscReal atol)
15934f43fa5SPieter Ghysels {
16034f43fa5SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
16134f43fa5SPieter Ghysels 
16234f43fa5SPieter Ghysels   PetscFunctionBegin;
163a36bf211SPieter Ghysels   PetscStackCall("STRUMPACK_set_HSS_abs_tol", STRUMPACK_set_HSS_abs_tol(*S,atol));
16434f43fa5SPieter Ghysels   PetscFunctionReturn(0);
16534f43fa5SPieter Ghysels }
166e5a36eccSStefano Zampini 
16734f43fa5SPieter Ghysels /*@
16834f43fa5SPieter Ghysels   MatSTRUMPACKSetHSSAbsTol - Set STRUMPACK absolute tolerance for HSS compression
169*147403d9SBarry Smith 
17034f43fa5SPieter Ghysels    Logically Collective on Mat
17134f43fa5SPieter Ghysels 
17234f43fa5SPieter Ghysels    Input Parameters:
17334f43fa5SPieter Ghysels +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
17434f43fa5SPieter Ghysels -  atol - absolute compression tolerance
17534f43fa5SPieter Ghysels 
17634f43fa5SPieter Ghysels   Options Database:
17734f43fa5SPieter Ghysels .   -mat_strumpack_hss_abs_tol <1e-8>         - Absolute compression tolerance (None)
17834f43fa5SPieter Ghysels 
17934f43fa5SPieter Ghysels    Level: beginner
18034f43fa5SPieter Ghysels 
18134f43fa5SPieter Ghysels    References:
182606c0280SSatish Balay .  * - STRUMPACK manual
18334f43fa5SPieter Ghysels 
18434f43fa5SPieter Ghysels .seealso: MatGetFactor()
18534f43fa5SPieter Ghysels @*/
18634f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSAbsTol(Mat F,PetscReal atol)
18734f43fa5SPieter Ghysels {
18834f43fa5SPieter Ghysels   PetscErrorCode ierr;
18934f43fa5SPieter Ghysels 
19034f43fa5SPieter Ghysels   PetscFunctionBegin;
19134f43fa5SPieter Ghysels   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
19234f43fa5SPieter Ghysels   PetscValidLogicalCollectiveReal(F,atol,2);
19334f43fa5SPieter Ghysels   ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSAbsTol_C",(Mat,PetscReal),(F,atol));CHKERRQ(ierr);
19434f43fa5SPieter Ghysels   PetscFunctionReturn(0);
19534f43fa5SPieter Ghysels }
19634f43fa5SPieter Ghysels 
19734f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSMaxRank_STRUMPACK(Mat F,PetscInt hssmaxrank)
19834f43fa5SPieter Ghysels {
19934f43fa5SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
20034f43fa5SPieter Ghysels 
20134f43fa5SPieter Ghysels   PetscFunctionBegin;
202a36bf211SPieter Ghysels   PetscStackCall("STRUMPACK_set_HSS_max_rank", STRUMPACK_set_HSS_max_rank(*S,hssmaxrank));
20334f43fa5SPieter Ghysels   PetscFunctionReturn(0);
20434f43fa5SPieter Ghysels }
205e5a36eccSStefano Zampini 
20634f43fa5SPieter Ghysels /*@
20734f43fa5SPieter Ghysels   MatSTRUMPACKSetHSSMaxRank - Set STRUMPACK maximum HSS rank
208*147403d9SBarry Smith 
20934f43fa5SPieter Ghysels    Logically Collective on Mat
21034f43fa5SPieter Ghysels 
21134f43fa5SPieter Ghysels    Input Parameters:
21234f43fa5SPieter Ghysels +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
21334f43fa5SPieter Ghysels -  hssmaxrank - maximum rank used in low-rank approximation
21434f43fa5SPieter Ghysels 
21534f43fa5SPieter Ghysels   Options Database:
21634f43fa5SPieter Ghysels .   -mat_strumpack_max_rank    - Maximum rank in HSS compression, when using pctype ilu (None)
21734f43fa5SPieter Ghysels 
21834f43fa5SPieter Ghysels    Level: beginner
21934f43fa5SPieter Ghysels 
22034f43fa5SPieter Ghysels    References:
221606c0280SSatish Balay .  * - STRUMPACK manual
22234f43fa5SPieter Ghysels 
22334f43fa5SPieter Ghysels .seealso: MatGetFactor()
22434f43fa5SPieter Ghysels @*/
22534f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSMaxRank(Mat F,PetscInt hssmaxrank)
22634f43fa5SPieter Ghysels {
22734f43fa5SPieter Ghysels   PetscErrorCode ierr;
22834f43fa5SPieter Ghysels 
22934f43fa5SPieter Ghysels   PetscFunctionBegin;
23034f43fa5SPieter Ghysels   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
23134f43fa5SPieter Ghysels   PetscValidLogicalCollectiveInt(F,hssmaxrank,2);
23234f43fa5SPieter Ghysels   ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSMaxRank_C",(Mat,PetscInt),(F,hssmaxrank));CHKERRQ(ierr);
23334f43fa5SPieter Ghysels   PetscFunctionReturn(0);
23434f43fa5SPieter Ghysels }
23534f43fa5SPieter Ghysels 
236a36bf211SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSLeafSize_STRUMPACK(Mat F,PetscInt leaf_size)
237a36bf211SPieter Ghysels {
238a36bf211SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
239a36bf211SPieter Ghysels 
240a36bf211SPieter Ghysels   PetscFunctionBegin;
241a36bf211SPieter Ghysels   PetscStackCall("STRUMPACK_set_HSS_leaf_size", STRUMPACK_set_HSS_leaf_size(*S,leaf_size));
242a36bf211SPieter Ghysels   PetscFunctionReturn(0);
243a36bf211SPieter Ghysels }
244e5a36eccSStefano Zampini 
245a36bf211SPieter Ghysels /*@
246a36bf211SPieter Ghysels   MatSTRUMPACKSetHSSLeafSize - Set STRUMPACK HSS leaf size
247*147403d9SBarry Smith 
248a36bf211SPieter Ghysels    Logically Collective on Mat
249a36bf211SPieter Ghysels 
250a36bf211SPieter Ghysels    Input Parameters:
251a36bf211SPieter Ghysels +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
252a36bf211SPieter Ghysels -  leaf_size - Size of diagonal blocks in HSS approximation
253a36bf211SPieter Ghysels 
254a36bf211SPieter Ghysels   Options Database:
255a36bf211SPieter Ghysels .   -mat_strumpack_leaf_size    - Size of diagonal blocks in HSS approximation, when using pctype ilu (None)
256a36bf211SPieter Ghysels 
257a36bf211SPieter Ghysels    Level: beginner
258a36bf211SPieter Ghysels 
259a36bf211SPieter Ghysels    References:
260606c0280SSatish Balay .  * - STRUMPACK manual
261a36bf211SPieter Ghysels 
262a36bf211SPieter Ghysels .seealso: MatGetFactor()
263a36bf211SPieter Ghysels @*/
264a36bf211SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSLeafSize(Mat F,PetscInt leaf_size)
265a36bf211SPieter Ghysels {
266a36bf211SPieter Ghysels   PetscErrorCode ierr;
267a36bf211SPieter Ghysels 
268a36bf211SPieter Ghysels   PetscFunctionBegin;
269a36bf211SPieter Ghysels   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
270a36bf211SPieter Ghysels   PetscValidLogicalCollectiveInt(F,leaf_size,2);
271a36bf211SPieter Ghysels   ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSLeafSize_C",(Mat,PetscInt),(F,leaf_size));CHKERRQ(ierr);
272a36bf211SPieter Ghysels   PetscFunctionReturn(0);
273a36bf211SPieter Ghysels }
274a36bf211SPieter Ghysels 
27534f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSMinSepSize_STRUMPACK(Mat F,PetscInt hssminsize)
27634f43fa5SPieter Ghysels {
27734f43fa5SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
27834f43fa5SPieter Ghysels 
27934f43fa5SPieter Ghysels   PetscFunctionBegin;
28034f43fa5SPieter Ghysels   PetscStackCall("STRUMPACK_set_HSS_min_sep_size", STRUMPACK_set_HSS_min_sep_size(*S,hssminsize));
28134f43fa5SPieter Ghysels   PetscFunctionReturn(0);
28234f43fa5SPieter Ghysels }
283e5a36eccSStefano Zampini 
284291d0ed5SPieter Ghysels /*@
285291d0ed5SPieter Ghysels   MatSTRUMPACKSetHSSMinSepSize - Set STRUMPACK minimum separator size for low-rank approximation
286*147403d9SBarry Smith 
287291d0ed5SPieter Ghysels    Logically Collective on Mat
288291d0ed5SPieter Ghysels 
289291d0ed5SPieter Ghysels    Input Parameters:
290291d0ed5SPieter Ghysels +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
291291d0ed5SPieter Ghysels -  hssminsize - minimum dense matrix size for low-rank approximation
292291d0ed5SPieter Ghysels 
293291d0ed5SPieter Ghysels   Options Database:
294*147403d9SBarry Smith .   -mat_strumpack_hss_min_sep_size <hssminsize> - set the minimum separator size
295291d0ed5SPieter Ghysels 
296291d0ed5SPieter Ghysels    Level: beginner
297291d0ed5SPieter Ghysels 
298291d0ed5SPieter Ghysels    References:
299606c0280SSatish Balay .  * - STRUMPACK manual
300291d0ed5SPieter Ghysels 
301291d0ed5SPieter Ghysels .seealso: MatGetFactor()
302291d0ed5SPieter Ghysels @*/
303291d0ed5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat F,PetscInt hssminsize)
304291d0ed5SPieter Ghysels {
305291d0ed5SPieter Ghysels   PetscErrorCode ierr;
306291d0ed5SPieter Ghysels 
307291d0ed5SPieter Ghysels   PetscFunctionBegin;
308291d0ed5SPieter Ghysels   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
309291d0ed5SPieter Ghysels   PetscValidLogicalCollectiveInt(F,hssminsize,2);
310291d0ed5SPieter Ghysels   ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSMinSepSize_C",(Mat,PetscInt),(F,hssminsize));CHKERRQ(ierr);
311291d0ed5SPieter Ghysels   PetscFunctionReturn(0);
312291d0ed5SPieter Ghysels }
313291d0ed5SPieter Ghysels 
314ad0c5e61SPieter Ghysels static PetscErrorCode MatSolve_STRUMPACK(Mat A,Vec b_mpi,Vec x)
3157d6ea485SPieter Ghysels {
3160d08a34dSPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)A->spptr;
3177d6ea485SPieter Ghysels   STRUMPACK_RETURN_CODE  sp_err;
3187d6ea485SPieter Ghysels   PetscErrorCode         ierr;
3197d6ea485SPieter Ghysels   const PetscScalar      *bptr;
3207d6ea485SPieter Ghysels   PetscScalar            *xptr;
3217d6ea485SPieter Ghysels 
3227d6ea485SPieter Ghysels   PetscFunctionBegin;
3237d6ea485SPieter Ghysels   ierr = VecGetArray(x,&xptr);CHKERRQ(ierr);
3247d6ea485SPieter Ghysels   ierr = VecGetArrayRead(b_mpi,&bptr);CHKERRQ(ierr);
3250d08a34dSPieter Ghysels 
3260d08a34dSPieter Ghysels   PetscStackCall("STRUMPACK_solve",sp_err = STRUMPACK_solve(*S,(PetscScalar*)bptr,xptr,0));
3270d08a34dSPieter Ghysels   switch (sp_err) {
3280d08a34dSPieter Ghysels   case STRUMPACK_SUCCESS: break;
3290d08a34dSPieter Ghysels   case STRUMPACK_MATRIX_NOT_SET:   { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); break; }
3300d08a34dSPieter Ghysels   case STRUMPACK_REORDERING_ERROR: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); break; }
3310d08a34dSPieter Ghysels   default:                           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: solve failed");
3327d6ea485SPieter Ghysels   }
3337d6ea485SPieter Ghysels   ierr = VecRestoreArray(x,&xptr);CHKERRQ(ierr);
3347d6ea485SPieter Ghysels   ierr = VecRestoreArrayRead(b_mpi,&bptr);CHKERRQ(ierr);
3357d6ea485SPieter Ghysels   PetscFunctionReturn(0);
3367d6ea485SPieter Ghysels }
3377d6ea485SPieter Ghysels 
338ad0c5e61SPieter Ghysels static PetscErrorCode MatMatSolve_STRUMPACK(Mat A,Mat B_mpi,Mat X)
3397d6ea485SPieter Ghysels {
3407d6ea485SPieter Ghysels   PetscErrorCode   ierr;
3417d6ea485SPieter Ghysels   PetscBool        flg;
3427d6ea485SPieter Ghysels 
3437d6ea485SPieter Ghysels   PetscFunctionBegin;
3447d6ea485SPieter Ghysels   ierr = PetscObjectTypeCompareAny((PetscObject)B_mpi,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
3452c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!flg,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
3467d6ea485SPieter Ghysels   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
3472c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!flg,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
3487d6ea485SPieter Ghysels   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_STRUMPACK() is not implemented yet");
3497d6ea485SPieter Ghysels   PetscFunctionReturn(0);
3507d6ea485SPieter Ghysels }
3517d6ea485SPieter Ghysels 
352860c79edSBarry Smith static PetscErrorCode MatView_Info_STRUMPACK(Mat A,PetscViewer viewer)
353ad0c5e61SPieter Ghysels {
354ad0c5e61SPieter Ghysels   PetscErrorCode  ierr;
355ad0c5e61SPieter Ghysels 
356ad0c5e61SPieter Ghysels   PetscFunctionBegin;
357ad0c5e61SPieter Ghysels   /* check if matrix is strumpack type */
358ad0c5e61SPieter Ghysels   if (A->ops->solve != MatSolve_STRUMPACK) PetscFunctionReturn(0);
359ad0c5e61SPieter Ghysels   ierr = PetscViewerASCIIPrintf(viewer,"STRUMPACK sparse solver!\n");CHKERRQ(ierr);
360ad0c5e61SPieter Ghysels   PetscFunctionReturn(0);
361ad0c5e61SPieter Ghysels }
362ad0c5e61SPieter Ghysels 
363ad0c5e61SPieter Ghysels static PetscErrorCode MatView_STRUMPACK(Mat A,PetscViewer viewer)
364ad0c5e61SPieter Ghysels {
365ad0c5e61SPieter Ghysels   PetscErrorCode    ierr;
366ad0c5e61SPieter Ghysels   PetscBool         iascii;
367ad0c5e61SPieter Ghysels   PetscViewerFormat format;
368ad0c5e61SPieter Ghysels 
369ad0c5e61SPieter Ghysels   PetscFunctionBegin;
370ad0c5e61SPieter Ghysels   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
371ad0c5e61SPieter Ghysels   if (iascii) {
372ad0c5e61SPieter Ghysels     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
373ad0c5e61SPieter Ghysels     if (format == PETSC_VIEWER_ASCII_INFO) {
374860c79edSBarry Smith       ierr = MatView_Info_STRUMPACK(A,viewer);CHKERRQ(ierr);
375ad0c5e61SPieter Ghysels     }
376ad0c5e61SPieter Ghysels   }
377ad0c5e61SPieter Ghysels   PetscFunctionReturn(0);
378ad0c5e61SPieter Ghysels }
3797d6ea485SPieter Ghysels 
380ad0c5e61SPieter Ghysels static PetscErrorCode MatLUFactorNumeric_STRUMPACK(Mat F,Mat A,const MatFactorInfo *info)
3817d6ea485SPieter Ghysels {
3820d08a34dSPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
3837d6ea485SPieter Ghysels   STRUMPACK_RETURN_CODE  sp_err;
3847d6ea485SPieter Ghysels   Mat_SeqAIJ             *A_d,*A_o;
3857d6ea485SPieter Ghysels   Mat_MPIAIJ             *mat;
3867d6ea485SPieter Ghysels   PetscErrorCode         ierr;
3870d08a34dSPieter Ghysels   PetscInt               M=A->rmap->N,m=A->rmap->n;
3887d6ea485SPieter Ghysels   PetscBool              flg;
3897d6ea485SPieter Ghysels 
3907d6ea485SPieter Ghysels   PetscFunctionBegin;
3917d6ea485SPieter Ghysels   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
3927d6ea485SPieter Ghysels   if (flg) { /* A is MATMPIAIJ */
3937d6ea485SPieter Ghysels     mat = (Mat_MPIAIJ*)A->data;
3947d6ea485SPieter Ghysels     A_d = (Mat_SeqAIJ*)(mat->A)->data;
3957d6ea485SPieter Ghysels     A_o = (Mat_SeqAIJ*)(mat->B)->data;
3960d08a34dSPieter Ghysels     PetscStackCall("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));
3977d6ea485SPieter Ghysels   } else { /* A is MATSEQAIJ */
3987d6ea485SPieter Ghysels     A_d = (Mat_SeqAIJ*)A->data;
3990d08a34dSPieter Ghysels     PetscStackCall("STRUMPACK_set_csr_matrix",STRUMPACK_set_csr_matrix(*S,&M,A_d->i,A_d->j,A_d->a,0));
4007d6ea485SPieter Ghysels   }
4017d6ea485SPieter Ghysels 
4027d6ea485SPieter Ghysels   /* Reorder and Factor the matrix. */
4037d6ea485SPieter Ghysels   /* TODO figure out how to avoid reorder if the matrix values changed, but the pattern remains the same. */
4040d08a34dSPieter Ghysels   PetscStackCall("STRUMPACK_reorder",sp_err = STRUMPACK_reorder(*S));
4050d08a34dSPieter Ghysels   PetscStackCall("STRUMPACK_factor",sp_err = STRUMPACK_factor(*S));
4060d08a34dSPieter Ghysels   switch (sp_err) {
4070d08a34dSPieter Ghysels   case STRUMPACK_SUCCESS: break;
4080d08a34dSPieter Ghysels   case STRUMPACK_MATRIX_NOT_SET:   { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); break; }
4090d08a34dSPieter Ghysels   case STRUMPACK_REORDERING_ERROR: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); break; }
4100d08a34dSPieter Ghysels   default:                           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: factorization failed");
4117d6ea485SPieter Ghysels   }
412cb250fa3SPieter Ghysels   F->assembled    = PETSC_TRUE;
413cb250fa3SPieter Ghysels   F->preallocated = PETSC_TRUE;
4147d6ea485SPieter Ghysels   PetscFunctionReturn(0);
4157d6ea485SPieter Ghysels }
4167d6ea485SPieter Ghysels 
417ad0c5e61SPieter Ghysels static PetscErrorCode MatLUFactorSymbolic_STRUMPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
4187d6ea485SPieter Ghysels {
4197d6ea485SPieter Ghysels   PetscFunctionBegin;
4207d6ea485SPieter Ghysels   F->ops->lufactornumeric = MatLUFactorNumeric_STRUMPACK;
4217d6ea485SPieter Ghysels   F->ops->solve           = MatSolve_STRUMPACK;
4227d6ea485SPieter Ghysels   F->ops->matsolve        = MatMatSolve_STRUMPACK;
4237d6ea485SPieter Ghysels   PetscFunctionReturn(0);
4247d6ea485SPieter Ghysels }
4257d6ea485SPieter Ghysels 
426ea799195SBarry Smith static PetscErrorCode MatFactorGetSolverType_aij_strumpack(Mat A,MatSolverType *type)
4277d6ea485SPieter Ghysels {
4287d6ea485SPieter Ghysels   PetscFunctionBegin;
4297d6ea485SPieter Ghysels   *type = MATSOLVERSTRUMPACK;
4307d6ea485SPieter Ghysels   PetscFunctionReturn(0);
4317d6ea485SPieter Ghysels }
4327d6ea485SPieter Ghysels 
433575704cbSPieter Ghysels /*MC
434575704cbSPieter Ghysels   MATSOLVERSSTRUMPACK = "strumpack" - A solver package providing a direct sparse solver (PCLU)
435575704cbSPieter Ghysels   and a preconditioner (PCILU) using low-rank compression via the external package STRUMPACK.
436575704cbSPieter Ghysels 
437575704cbSPieter Ghysels   Consult the STRUMPACK-sparse manual for more info.
438575704cbSPieter Ghysels 
439575704cbSPieter Ghysels   Use
440575704cbSPieter Ghysels      ./configure --download-strumpack
441575704cbSPieter Ghysels   to have PETSc installed with STRUMPACK
442575704cbSPieter Ghysels 
443575704cbSPieter Ghysels   Use
4443ca39a21SBarry Smith     -pc_type lu -pc_factor_mat_solver_type strumpack
445575704cbSPieter Ghysels   to use this as an exact (direct) solver, use
4463ca39a21SBarry Smith     -pc_type ilu -pc_factor_mat_solver_type strumpack
447575704cbSPieter Ghysels   to enable low-rank compression (i.e, use as a preconditioner).
448575704cbSPieter Ghysels 
449575704cbSPieter Ghysels   Works with AIJ matrices
450575704cbSPieter Ghysels 
451575704cbSPieter Ghysels   Options Database Keys:
45234f43fa5SPieter Ghysels + -mat_strumpack_verbose
45334f43fa5SPieter Ghysels . -mat_strumpack_hss_rel_tol <1e-2>         - Relative compression tolerance (None)
45434f43fa5SPieter Ghysels . -mat_strumpack_hss_abs_tol <1e-8>         - Absolute compression tolerance (None)
455575704cbSPieter Ghysels . -mat_strumpack_colperm <TRUE>             - Permute matrix to make diagonal nonzeros (None)
45634f43fa5SPieter Ghysels . -mat_strumpack_hss_min_sep_size <256>     - Minimum size of separator for HSS compression (None)
45734f43fa5SPieter Ghysels . -mat_strumpack_max_rank                   - Maximum rank in HSS compression, when using pctype ilu (None)
458a36bf211SPieter Ghysels . -mat_strumpack_leaf_size                  - Size of diagonal blocks in HSS approximation, when using pctype ilu (None)
45934f43fa5SPieter Ghysels . -mat_strumpack_reordering <METIS>         - Sparsity reducing matrix reordering (choose one of) NATURAL METIS PARMETIS SCOTCH PTSCOTCH RCM (None)
460a2b725a8SWilliam Gropp - -mat_strumpack_iterative_solver <DIRECT>  - Select iterative solver from STRUMPACK (choose one of) AUTO DIRECT REFINE PREC_GMRES GMRES PREC_BICGSTAB BICGSTAB (None)
461575704cbSPieter Ghysels 
462575704cbSPieter Ghysels  Level: beginner
463575704cbSPieter Ghysels 
4643ca39a21SBarry Smith .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, PCFactorSetMatSolverType(), MatSolverType
465575704cbSPieter Ghysels M*/
466ad0c5e61SPieter Ghysels static PetscErrorCode MatGetFactor_aij_strumpack(Mat A,MatFactorType ftype,Mat *F)
4677d6ea485SPieter Ghysels {
4687d6ea485SPieter Ghysels   Mat                           B;
4697d6ea485SPieter Ghysels   PetscErrorCode                ierr;
4707d6ea485SPieter Ghysels   PetscInt                      M=A->rmap->N,N=A->cmap->N;
471575704cbSPieter Ghysels   PetscBool                     verb,flg,set;
47234f43fa5SPieter Ghysels   PetscReal                     ctol;
473a36bf211SPieter Ghysels   PetscInt                      hssminsize,max_rank,leaf_size;
47434f43fa5SPieter Ghysels   STRUMPACK_SparseSolver        *S;
47534f43fa5SPieter Ghysels   STRUMPACK_INTERFACE           iface;
47634f43fa5SPieter Ghysels   STRUMPACK_REORDERING_STRATEGY ndcurrent,ndvalue;
47734f43fa5SPieter Ghysels   STRUMPACK_KRYLOV_SOLVER       itcurrent,itsolver;
47834f43fa5SPieter Ghysels   const STRUMPACK_PRECISION     table[2][2][2] =
47934f43fa5SPieter Ghysels     {{{STRUMPACK_FLOATCOMPLEX_64, STRUMPACK_DOUBLECOMPLEX_64},
4807d6ea485SPieter Ghysels       {STRUMPACK_FLOAT_64,        STRUMPACK_DOUBLE_64}},
4817d6ea485SPieter Ghysels      {{STRUMPACK_FLOATCOMPLEX,    STRUMPACK_DOUBLECOMPLEX},
4827d6ea485SPieter Ghysels       {STRUMPACK_FLOAT,           STRUMPACK_DOUBLE}}};
4834ac6704cSBarry Smith   const STRUMPACK_PRECISION     prec = table[(sizeof(PetscInt)==8)?0:1][(PETSC_SCALAR==PETSC_COMPLEX)?0:1][(PETSC_REAL==PETSC_FLOAT)?0:1];
4844ac6704cSBarry Smith   const char *const             STRUMPACKNDTypes[] = {"NATURAL","METIS","PARMETIS","SCOTCH","PTSCOTCH","RCM","STRUMPACKNDTypes","",0};
4854ac6704cSBarry Smith   const char *const             SolverTypes[] = {"AUTO","NONE","REFINE","PREC_GMRES","GMRES","PREC_BICGSTAB","BICGSTAB","SolverTypes","",0};
4867d6ea485SPieter Ghysels 
4877d6ea485SPieter Ghysels   PetscFunctionBegin;
4887d6ea485SPieter Ghysels   /* Create the factorization matrix */
4897d6ea485SPieter Ghysels   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
4907d6ea485SPieter Ghysels   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,M,N);CHKERRQ(ierr);
4917d6ea485SPieter Ghysels   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
4921e1ea65dSPierre Jolivet   ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr);
4937d6ea485SPieter Ghysels   ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr);
49466e17bc3SBarry Smith   B->trivialsymbolic = PETSC_TRUE;
495575704cbSPieter Ghysels   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) {
4967d6ea485SPieter Ghysels     B->ops->lufactorsymbolic  = MatLUFactorSymbolic_STRUMPACK;
497575704cbSPieter Ghysels     B->ops->ilufactorsymbolic = MatLUFactorSymbolic_STRUMPACK;
498575704cbSPieter Ghysels   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
4997d6ea485SPieter Ghysels   B->ops->view        = MatView_STRUMPACK;
5007d6ea485SPieter Ghysels   B->ops->destroy     = MatDestroy_STRUMPACK;
5017d6ea485SPieter Ghysels   B->ops->getdiagonal = MatGetDiagonal_STRUMPACK;
5023ca39a21SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_aij_strumpack);CHKERRQ(ierr);
50334f43fa5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetReordering_C",MatSTRUMPACKSetReordering_STRUMPACK);CHKERRQ(ierr);
50434f43fa5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetColPerm_C",MatSTRUMPACKSetColPerm_STRUMPACK);CHKERRQ(ierr);
50534f43fa5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSRelTol_C",MatSTRUMPACKSetHSSRelTol_STRUMPACK);CHKERRQ(ierr);
50634f43fa5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSAbsTol_C",MatSTRUMPACKSetHSSAbsTol_STRUMPACK);CHKERRQ(ierr);
507a36bf211SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSMaxRank_C",MatSTRUMPACKSetHSSMaxRank_STRUMPACK);CHKERRQ(ierr);
508a36bf211SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSLeafSize_C",MatSTRUMPACKSetHSSLeafSize_STRUMPACK);CHKERRQ(ierr);
509291d0ed5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSMinSepSize_C",MatSTRUMPACKSetHSSMinSepSize_STRUMPACK);CHKERRQ(ierr);
510575704cbSPieter Ghysels   B->factortype = ftype;
5110d08a34dSPieter Ghysels   ierr     = PetscNewLog(B,&S);CHKERRQ(ierr);
5120d08a34dSPieter Ghysels   B->spptr = S;
5130d08a34dSPieter Ghysels 
5141e1ea65dSPierre Jolivet   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr);
5150d08a34dSPieter Ghysels   iface = flg ? STRUMPACK_MT : STRUMPACK_MPI_DIST;
5167d6ea485SPieter Ghysels 
5177d6ea485SPieter Ghysels   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"STRUMPACK Options","Mat");CHKERRQ(ierr);
5187d6ea485SPieter Ghysels 
519575704cbSPieter Ghysels   verb = PetscLogPrintInfo ? PETSC_TRUE : PETSC_FALSE;
5207d6ea485SPieter Ghysels   ierr = PetscOptionsBool("-mat_strumpack_verbose","Print STRUMPACK information","None",verb,&verb,NULL);CHKERRQ(ierr);
5217d6ea485SPieter Ghysels 
52234f43fa5SPieter Ghysels   PetscStackCall("STRUMPACK_init",STRUMPACK_init(S,PetscObjectComm((PetscObject)A),prec,iface,0,NULL,verb));
52355c022e5SPieter Ghysels 
524a36bf211SPieter Ghysels   PetscStackCall("STRUMPACK_HSS_rel_tol",ctol = (PetscReal)STRUMPACK_HSS_rel_tol(*S));
52534f43fa5SPieter Ghysels   ierr = PetscOptionsReal("-mat_strumpack_hss_rel_tol","Relative compression tolerance","None",ctol,&ctol,&set);CHKERRQ(ierr);
526a36bf211SPieter Ghysels   if (set) PetscStackCall("STRUMPACK_set_HSS_rel_tol",STRUMPACK_set_HSS_rel_tol(*S,(double)ctol));
5277d6ea485SPieter Ghysels 
528a36bf211SPieter Ghysels   PetscStackCall("STRUMPACK_HSS_abs_tol",ctol = (PetscReal)STRUMPACK_HSS_abs_tol(*S));
52934f43fa5SPieter Ghysels   ierr = PetscOptionsReal("-mat_strumpack_hss_abs_tol","Absolute compression tolerance","None",ctol,&ctol,&set);CHKERRQ(ierr);
530a36bf211SPieter Ghysels   if (set) PetscStackCall("STRUMPACK_set_HSS_abs_tol",STRUMPACK_set_HSS_abs_tol(*S,(double)ctol));
531575704cbSPieter Ghysels 
532291d0ed5SPieter Ghysels   PetscStackCall("STRUMPACK_mc64job",flg = (STRUMPACK_mc64job(*S) == 0) ? PETSC_FALSE : PETSC_TRUE);
533575704cbSPieter Ghysels   ierr = PetscOptionsBool("-mat_strumpack_colperm","Find a col perm to get nonzero diagonal","None",flg,&flg,&set);CHKERRQ(ierr);
5340d08a34dSPieter Ghysels   if (set) PetscStackCall("STRUMPACK_set_mc64job",STRUMPACK_set_mc64job(*S,flg ? 5 : 0));
535575704cbSPieter Ghysels 
536291d0ed5SPieter Ghysels   PetscStackCall("STRUMPACK_HSS_min_sep_size",hssminsize = (PetscInt)STRUMPACK_HSS_min_sep_size(*S));
537291d0ed5SPieter Ghysels   ierr = PetscOptionsInt("-mat_strumpack_hss_min_sep_size","Minimum size of separator for HSS compression","None",hssminsize,&hssminsize,&set);CHKERRQ(ierr);
538291d0ed5SPieter Ghysels   if (set) PetscStackCall("STRUMPACK_set_HSS_min_sep_size",STRUMPACK_set_HSS_min_sep_size(*S,(int)hssminsize));
539575704cbSPieter Ghysels 
540a36bf211SPieter Ghysels   PetscStackCall("STRUMPACK_HSS_max_rank",max_rank = (PetscInt)STRUMPACK_HSS_max_rank(*S));
54134f43fa5SPieter Ghysels   ierr = PetscOptionsInt("-mat_strumpack_max_rank","Maximum rank in HSS approximation","None",max_rank,&max_rank,&set);CHKERRQ(ierr);
542a36bf211SPieter Ghysels   if (set) PetscStackCall("STRUMPACK_set_HSS_max_rank",STRUMPACK_set_HSS_max_rank(*S,(int)max_rank));
543a36bf211SPieter Ghysels 
544a36bf211SPieter Ghysels   PetscStackCall("STRUMPACK_HSS_leaf_size",leaf_size = (PetscInt)STRUMPACK_HSS_leaf_size(*S));
545a36bf211SPieter Ghysels   ierr = PetscOptionsInt("-mat_strumpack_leaf_size","Size of diagonal blocks in HSS approximation","None",leaf_size,&leaf_size,&set);CHKERRQ(ierr);
546a36bf211SPieter Ghysels   if (set) PetscStackCall("STRUMPACK_set_HSS_leaf_size",STRUMPACK_set_HSS_leaf_size(*S,(int)leaf_size));
54734f43fa5SPieter Ghysels 
54834f43fa5SPieter Ghysels   PetscStackCall("STRUMPACK_reordering_method",ndcurrent = STRUMPACK_reordering_method(*S));
5492f613bf5SBarry Smith   ierr = PetscOptionsEnum("-mat_strumpack_reordering","Sparsity reducing matrix reordering","None",STRUMPACKNDTypes,(PetscEnum)ndcurrent,(PetscEnum*)&ndvalue,&set);CHKERRQ(ierr);
55034f43fa5SPieter Ghysels   if (set) PetscStackCall("STRUMPACK_set_reordering_method",STRUMPACK_set_reordering_method(*S,ndvalue));
551575704cbSPieter Ghysels 
55288382b05SPieter Ghysels   if (ftype == MAT_FACTOR_ILU) {
55388382b05SPieter Ghysels     /* When enabling HSS compression, the STRUMPACK solver becomes an incomplete              */
55488382b05SPieter Ghysels     /* (or approximate) LU factorization.                                                     */
555291d0ed5SPieter Ghysels     PetscStackCall("STRUMPACK_enable_HSS",STRUMPACK_enable_HSS(*S));
55688382b05SPieter Ghysels   }
557575704cbSPieter Ghysels 
55834f43fa5SPieter Ghysels   /* Disable the outer iterative solver from STRUMPACK.                                       */
55934f43fa5SPieter Ghysels   /* When STRUMPACK is used as a direct solver, it will by default do iterative refinement.   */
56034f43fa5SPieter Ghysels   /* When STRUMPACK is used as an approximate factorization preconditioner (by enabling       */
56134f43fa5SPieter Ghysels   /* low-rank compression), it will use it's own preconditioned GMRES. Here we can disable    */
56234f43fa5SPieter Ghysels   /* the outer iterative solver, as PETSc uses STRUMPACK from within a KSP.                   */
56334f43fa5SPieter Ghysels   PetscStackCall("STRUMPACK_set_Krylov_solver", STRUMPACK_set_Krylov_solver(*S, STRUMPACK_DIRECT));
56434f43fa5SPieter Ghysels 
56534f43fa5SPieter Ghysels   PetscStackCall("STRUMPACK_Krylov_solver",itcurrent = STRUMPACK_Krylov_solver(*S));
5662f613bf5SBarry Smith   ierr = PetscOptionsEnum("-mat_strumpack_iterative_solver","Select iterative solver from STRUMPACK","None",SolverTypes,(PetscEnum)itcurrent,(PetscEnum*)&itsolver,&set);CHKERRQ(ierr);
56734f43fa5SPieter Ghysels   if (set) PetscStackCall("STRUMPACK_set_Krylov_solver",STRUMPACK_set_Krylov_solver(*S,itsolver));
56834f43fa5SPieter Ghysels 
56934f43fa5SPieter Ghysels   PetscOptionsEnd();
57055c022e5SPieter Ghysels 
5717d6ea485SPieter Ghysels   *F = B;
5727d6ea485SPieter Ghysels   PetscFunctionReturn(0);
5737d6ea485SPieter Ghysels }
5747d6ea485SPieter Ghysels 
5753ca39a21SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_STRUMPACK(void)
5767d6ea485SPieter Ghysels {
5777d6ea485SPieter Ghysels   PetscErrorCode ierr;
5787d6ea485SPieter Ghysels 
5797d6ea485SPieter Ghysels   PetscFunctionBegin;
5803ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr);
5813ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr);
5823ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_ILU,MatGetFactor_aij_strumpack);CHKERRQ(ierr);
5833ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_ILU,MatGetFactor_aij_strumpack);CHKERRQ(ierr);
5847d6ea485SPieter Ghysels   PetscFunctionReturn(0);
5857d6ea485SPieter Ghysels }
586