xref: /petsc/src/mat/impls/aij/mpi/strumpack/strumpack.c (revision 2f613bf53f46f9356e00a2ca2bd69453be72fc31)
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:
6534f43fa5SPieter Ghysels .      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
91575704cbSPieter Ghysels    Logically Collective on Mat
92575704cbSPieter Ghysels 
93575704cbSPieter Ghysels    Input Parameters:
94575704cbSPieter Ghysels +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
9534f43fa5SPieter Ghysels -  cperm - whether or not to permute (internally) the columns of the matrix
96575704cbSPieter Ghysels 
97575704cbSPieter Ghysels   Options Database:
9834f43fa5SPieter Ghysels .   -mat_strumpack_colperm <cperm>
99575704cbSPieter Ghysels 
100575704cbSPieter Ghysels    Level: beginner
101575704cbSPieter Ghysels 
102575704cbSPieter Ghysels    References:
103575704cbSPieter Ghysels .      STRUMPACK manual
104575704cbSPieter Ghysels 
105575704cbSPieter Ghysels .seealso: MatGetFactor()
106575704cbSPieter Ghysels @*/
10734f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetColPerm(Mat F,PetscBool cperm)
108575704cbSPieter Ghysels {
109575704cbSPieter Ghysels   PetscErrorCode ierr;
110575704cbSPieter Ghysels 
111575704cbSPieter Ghysels   PetscFunctionBegin;
112575704cbSPieter Ghysels   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
11334f43fa5SPieter Ghysels   PetscValidLogicalCollectiveBool(F,cperm,2);
11434f43fa5SPieter Ghysels   ierr = PetscTryMethod(F,"MatSTRUMPACKSetColPerm_C",(Mat,PetscBool),(F,cperm));CHKERRQ(ierr);
115575704cbSPieter Ghysels   PetscFunctionReturn(0);
116575704cbSPieter Ghysels }
117575704cbSPieter Ghysels 
11834f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSRelTol_STRUMPACK(Mat F,PetscReal rtol)
11934f43fa5SPieter Ghysels {
12034f43fa5SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
12134f43fa5SPieter Ghysels 
12234f43fa5SPieter Ghysels   PetscFunctionBegin;
123a36bf211SPieter Ghysels   PetscStackCall("STRUMPACK_set_HSS_rel_tol", STRUMPACK_set_HSS_rel_tol(*S,rtol));
12434f43fa5SPieter Ghysels   PetscFunctionReturn(0);
12534f43fa5SPieter Ghysels }
126e5a36eccSStefano Zampini 
12734f43fa5SPieter Ghysels /*@
12834f43fa5SPieter Ghysels   MatSTRUMPACKSetHSSRelTol - Set STRUMPACK relative tolerance for HSS compression
12934f43fa5SPieter Ghysels    Logically Collective on Mat
13034f43fa5SPieter Ghysels 
13134f43fa5SPieter Ghysels    Input Parameters:
13234f43fa5SPieter Ghysels +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
13334f43fa5SPieter Ghysels -  rtol - relative compression tolerance
13434f43fa5SPieter Ghysels 
13534f43fa5SPieter Ghysels   Options Database:
13634f43fa5SPieter Ghysels .   -mat_strumpack_hss_rel_tol <1e-2>         - Relative compression tolerance (None)
13734f43fa5SPieter Ghysels 
13834f43fa5SPieter Ghysels    Level: beginner
13934f43fa5SPieter Ghysels 
14034f43fa5SPieter Ghysels    References:
14134f43fa5SPieter Ghysels .      STRUMPACK manual
14234f43fa5SPieter Ghysels 
14334f43fa5SPieter Ghysels .seealso: MatGetFactor()
14434f43fa5SPieter Ghysels @*/
14534f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSRelTol(Mat F,PetscReal rtol)
14634f43fa5SPieter Ghysels {
14734f43fa5SPieter Ghysels   PetscErrorCode ierr;
14834f43fa5SPieter Ghysels 
14934f43fa5SPieter Ghysels   PetscFunctionBegin;
15034f43fa5SPieter Ghysels   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
15134f43fa5SPieter Ghysels   PetscValidLogicalCollectiveReal(F,rtol,2);
15234f43fa5SPieter Ghysels   ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSRelTol_C",(Mat,PetscReal),(F,rtol));CHKERRQ(ierr);
15334f43fa5SPieter Ghysels   PetscFunctionReturn(0);
15434f43fa5SPieter Ghysels }
15534f43fa5SPieter Ghysels 
15634f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSAbsTol_STRUMPACK(Mat F,PetscReal atol)
15734f43fa5SPieter Ghysels {
15834f43fa5SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
15934f43fa5SPieter Ghysels 
16034f43fa5SPieter Ghysels   PetscFunctionBegin;
161a36bf211SPieter Ghysels   PetscStackCall("STRUMPACK_set_HSS_abs_tol", STRUMPACK_set_HSS_abs_tol(*S,atol));
16234f43fa5SPieter Ghysels   PetscFunctionReturn(0);
16334f43fa5SPieter Ghysels }
164e5a36eccSStefano Zampini 
16534f43fa5SPieter Ghysels /*@
16634f43fa5SPieter Ghysels   MatSTRUMPACKSetHSSAbsTol - Set STRUMPACK absolute tolerance for HSS compression
16734f43fa5SPieter Ghysels    Logically Collective on Mat
16834f43fa5SPieter Ghysels 
16934f43fa5SPieter Ghysels    Input Parameters:
17034f43fa5SPieter Ghysels +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
17134f43fa5SPieter Ghysels -  atol - absolute compression tolerance
17234f43fa5SPieter Ghysels 
17334f43fa5SPieter Ghysels   Options Database:
17434f43fa5SPieter Ghysels .   -mat_strumpack_hss_abs_tol <1e-8>         - Absolute compression tolerance (None)
17534f43fa5SPieter Ghysels 
17634f43fa5SPieter Ghysels    Level: beginner
17734f43fa5SPieter Ghysels 
17834f43fa5SPieter Ghysels    References:
17934f43fa5SPieter Ghysels .      STRUMPACK manual
18034f43fa5SPieter Ghysels 
18134f43fa5SPieter Ghysels .seealso: MatGetFactor()
18234f43fa5SPieter Ghysels @*/
18334f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSAbsTol(Mat F,PetscReal atol)
18434f43fa5SPieter Ghysels {
18534f43fa5SPieter Ghysels   PetscErrorCode ierr;
18634f43fa5SPieter Ghysels 
18734f43fa5SPieter Ghysels   PetscFunctionBegin;
18834f43fa5SPieter Ghysels   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
18934f43fa5SPieter Ghysels   PetscValidLogicalCollectiveReal(F,atol,2);
19034f43fa5SPieter Ghysels   ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSAbsTol_C",(Mat,PetscReal),(F,atol));CHKERRQ(ierr);
19134f43fa5SPieter Ghysels   PetscFunctionReturn(0);
19234f43fa5SPieter Ghysels }
19334f43fa5SPieter Ghysels 
19434f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSMaxRank_STRUMPACK(Mat F,PetscInt hssmaxrank)
19534f43fa5SPieter Ghysels {
19634f43fa5SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
19734f43fa5SPieter Ghysels 
19834f43fa5SPieter Ghysels   PetscFunctionBegin;
199a36bf211SPieter Ghysels   PetscStackCall("STRUMPACK_set_HSS_max_rank", STRUMPACK_set_HSS_max_rank(*S,hssmaxrank));
20034f43fa5SPieter Ghysels   PetscFunctionReturn(0);
20134f43fa5SPieter Ghysels }
202e5a36eccSStefano Zampini 
20334f43fa5SPieter Ghysels /*@
20434f43fa5SPieter Ghysels   MatSTRUMPACKSetHSSMaxRank - Set STRUMPACK maximum HSS rank
20534f43fa5SPieter Ghysels    Logically Collective on Mat
20634f43fa5SPieter Ghysels 
20734f43fa5SPieter Ghysels    Input Parameters:
20834f43fa5SPieter Ghysels +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
20934f43fa5SPieter Ghysels -  hssmaxrank - maximum rank used in low-rank approximation
21034f43fa5SPieter Ghysels 
21134f43fa5SPieter Ghysels   Options Database:
21234f43fa5SPieter Ghysels .   -mat_strumpack_max_rank    - Maximum rank in HSS compression, when using pctype ilu (None)
21334f43fa5SPieter Ghysels 
21434f43fa5SPieter Ghysels    Level: beginner
21534f43fa5SPieter Ghysels 
21634f43fa5SPieter Ghysels    References:
21734f43fa5SPieter Ghysels .      STRUMPACK manual
21834f43fa5SPieter Ghysels 
21934f43fa5SPieter Ghysels .seealso: MatGetFactor()
22034f43fa5SPieter Ghysels @*/
22134f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSMaxRank(Mat F,PetscInt hssmaxrank)
22234f43fa5SPieter Ghysels {
22334f43fa5SPieter Ghysels   PetscErrorCode ierr;
22434f43fa5SPieter Ghysels 
22534f43fa5SPieter Ghysels   PetscFunctionBegin;
22634f43fa5SPieter Ghysels   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
22734f43fa5SPieter Ghysels   PetscValidLogicalCollectiveInt(F,hssmaxrank,2);
22834f43fa5SPieter Ghysels   ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSMaxRank_C",(Mat,PetscInt),(F,hssmaxrank));CHKERRQ(ierr);
22934f43fa5SPieter Ghysels   PetscFunctionReturn(0);
23034f43fa5SPieter Ghysels }
23134f43fa5SPieter Ghysels 
232a36bf211SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSLeafSize_STRUMPACK(Mat F,PetscInt leaf_size)
233a36bf211SPieter Ghysels {
234a36bf211SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
235a36bf211SPieter Ghysels 
236a36bf211SPieter Ghysels   PetscFunctionBegin;
237a36bf211SPieter Ghysels   PetscStackCall("STRUMPACK_set_HSS_leaf_size", STRUMPACK_set_HSS_leaf_size(*S,leaf_size));
238a36bf211SPieter Ghysels   PetscFunctionReturn(0);
239a36bf211SPieter Ghysels }
240e5a36eccSStefano Zampini 
241a36bf211SPieter Ghysels /*@
242a36bf211SPieter Ghysels   MatSTRUMPACKSetHSSLeafSize - Set STRUMPACK HSS leaf size
243a36bf211SPieter Ghysels    Logically Collective on Mat
244a36bf211SPieter Ghysels 
245a36bf211SPieter Ghysels    Input Parameters:
246a36bf211SPieter Ghysels +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
247a36bf211SPieter Ghysels -  leaf_size - Size of diagonal blocks in HSS approximation
248a36bf211SPieter Ghysels 
249a36bf211SPieter Ghysels   Options Database:
250a36bf211SPieter Ghysels .   -mat_strumpack_leaf_size    - Size of diagonal blocks in HSS approximation, when using pctype ilu (None)
251a36bf211SPieter Ghysels 
252a36bf211SPieter Ghysels    Level: beginner
253a36bf211SPieter Ghysels 
254a36bf211SPieter Ghysels    References:
255a36bf211SPieter Ghysels .      STRUMPACK manual
256a36bf211SPieter Ghysels 
257a36bf211SPieter Ghysels .seealso: MatGetFactor()
258a36bf211SPieter Ghysels @*/
259a36bf211SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSLeafSize(Mat F,PetscInt leaf_size)
260a36bf211SPieter Ghysels {
261a36bf211SPieter Ghysels   PetscErrorCode ierr;
262a36bf211SPieter Ghysels 
263a36bf211SPieter Ghysels   PetscFunctionBegin;
264a36bf211SPieter Ghysels   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
265a36bf211SPieter Ghysels   PetscValidLogicalCollectiveInt(F,leaf_size,2);
266a36bf211SPieter Ghysels   ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSLeafSize_C",(Mat,PetscInt),(F,leaf_size));CHKERRQ(ierr);
267a36bf211SPieter Ghysels   PetscFunctionReturn(0);
268a36bf211SPieter Ghysels }
269a36bf211SPieter Ghysels 
27034f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSMinSepSize_STRUMPACK(Mat F,PetscInt hssminsize)
27134f43fa5SPieter Ghysels {
27234f43fa5SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
27334f43fa5SPieter Ghysels 
27434f43fa5SPieter Ghysels   PetscFunctionBegin;
27534f43fa5SPieter Ghysels   PetscStackCall("STRUMPACK_set_HSS_min_sep_size", STRUMPACK_set_HSS_min_sep_size(*S,hssminsize));
27634f43fa5SPieter Ghysels   PetscFunctionReturn(0);
27734f43fa5SPieter Ghysels }
278e5a36eccSStefano Zampini 
279291d0ed5SPieter Ghysels /*@
280291d0ed5SPieter Ghysels   MatSTRUMPACKSetHSSMinSepSize - Set STRUMPACK minimum separator size for low-rank approximation
281291d0ed5SPieter Ghysels    Logically Collective on Mat
282291d0ed5SPieter Ghysels 
283291d0ed5SPieter Ghysels    Input Parameters:
284291d0ed5SPieter Ghysels +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
285291d0ed5SPieter Ghysels -  hssminsize - minimum dense matrix size for low-rank approximation
286291d0ed5SPieter Ghysels 
287291d0ed5SPieter Ghysels   Options Database:
288291d0ed5SPieter Ghysels .   -mat_strumpack_hss_min_sep_size <hssminsize>
289291d0ed5SPieter Ghysels 
290291d0ed5SPieter Ghysels    Level: beginner
291291d0ed5SPieter Ghysels 
292291d0ed5SPieter Ghysels    References:
293291d0ed5SPieter Ghysels .      STRUMPACK manual
294291d0ed5SPieter Ghysels 
295291d0ed5SPieter Ghysels .seealso: MatGetFactor()
296291d0ed5SPieter Ghysels @*/
297291d0ed5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat F,PetscInt hssminsize)
298291d0ed5SPieter Ghysels {
299291d0ed5SPieter Ghysels   PetscErrorCode ierr;
300291d0ed5SPieter Ghysels 
301291d0ed5SPieter Ghysels   PetscFunctionBegin;
302291d0ed5SPieter Ghysels   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
303291d0ed5SPieter Ghysels   PetscValidLogicalCollectiveInt(F,hssminsize,2);
304291d0ed5SPieter Ghysels   ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSMinSepSize_C",(Mat,PetscInt),(F,hssminsize));CHKERRQ(ierr);
305291d0ed5SPieter Ghysels   PetscFunctionReturn(0);
306291d0ed5SPieter Ghysels }
307291d0ed5SPieter Ghysels 
308ad0c5e61SPieter Ghysels static PetscErrorCode MatSolve_STRUMPACK(Mat A,Vec b_mpi,Vec x)
3097d6ea485SPieter Ghysels {
3100d08a34dSPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)A->spptr;
3117d6ea485SPieter Ghysels   STRUMPACK_RETURN_CODE  sp_err;
3127d6ea485SPieter Ghysels   PetscErrorCode         ierr;
3137d6ea485SPieter Ghysels   const PetscScalar      *bptr;
3147d6ea485SPieter Ghysels   PetscScalar            *xptr;
3157d6ea485SPieter Ghysels 
3167d6ea485SPieter Ghysels   PetscFunctionBegin;
3177d6ea485SPieter Ghysels   ierr = VecGetArray(x,&xptr);CHKERRQ(ierr);
3187d6ea485SPieter Ghysels   ierr = VecGetArrayRead(b_mpi,&bptr);CHKERRQ(ierr);
3190d08a34dSPieter Ghysels 
3200d08a34dSPieter Ghysels   PetscStackCall("STRUMPACK_solve",sp_err = STRUMPACK_solve(*S,(PetscScalar*)bptr,xptr,0));
3210d08a34dSPieter Ghysels   switch (sp_err) {
3220d08a34dSPieter Ghysels   case STRUMPACK_SUCCESS: break;
3230d08a34dSPieter Ghysels   case STRUMPACK_MATRIX_NOT_SET:   { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); break; }
3240d08a34dSPieter Ghysels   case STRUMPACK_REORDERING_ERROR: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); break; }
3250d08a34dSPieter Ghysels   default:                           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: solve failed");
3267d6ea485SPieter Ghysels   }
3277d6ea485SPieter Ghysels   ierr = VecRestoreArray(x,&xptr);CHKERRQ(ierr);
3287d6ea485SPieter Ghysels   ierr = VecRestoreArrayRead(b_mpi,&bptr);CHKERRQ(ierr);
3297d6ea485SPieter Ghysels   PetscFunctionReturn(0);
3307d6ea485SPieter Ghysels }
3317d6ea485SPieter Ghysels 
332ad0c5e61SPieter Ghysels static PetscErrorCode MatMatSolve_STRUMPACK(Mat A,Mat B_mpi,Mat X)
3337d6ea485SPieter Ghysels {
3347d6ea485SPieter Ghysels   PetscErrorCode   ierr;
3357d6ea485SPieter Ghysels   PetscBool        flg;
3367d6ea485SPieter Ghysels 
3377d6ea485SPieter Ghysels   PetscFunctionBegin;
3387d6ea485SPieter Ghysels   ierr = PetscObjectTypeCompareAny((PetscObject)B_mpi,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
3397d6ea485SPieter Ghysels   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
3407d6ea485SPieter Ghysels   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
3417d6ea485SPieter Ghysels   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
3427d6ea485SPieter Ghysels   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_STRUMPACK() is not implemented yet");
3437d6ea485SPieter Ghysels   PetscFunctionReturn(0);
3447d6ea485SPieter Ghysels }
3457d6ea485SPieter Ghysels 
346860c79edSBarry Smith static PetscErrorCode MatView_Info_STRUMPACK(Mat A,PetscViewer viewer)
347ad0c5e61SPieter Ghysels {
348ad0c5e61SPieter Ghysels   PetscErrorCode  ierr;
349ad0c5e61SPieter Ghysels 
350ad0c5e61SPieter Ghysels   PetscFunctionBegin;
351ad0c5e61SPieter Ghysels   /* check if matrix is strumpack type */
352ad0c5e61SPieter Ghysels   if (A->ops->solve != MatSolve_STRUMPACK) PetscFunctionReturn(0);
353ad0c5e61SPieter Ghysels   ierr = PetscViewerASCIIPrintf(viewer,"STRUMPACK sparse solver!\n");CHKERRQ(ierr);
354ad0c5e61SPieter Ghysels   PetscFunctionReturn(0);
355ad0c5e61SPieter Ghysels }
356ad0c5e61SPieter Ghysels 
357ad0c5e61SPieter Ghysels static PetscErrorCode MatView_STRUMPACK(Mat A,PetscViewer viewer)
358ad0c5e61SPieter Ghysels {
359ad0c5e61SPieter Ghysels   PetscErrorCode    ierr;
360ad0c5e61SPieter Ghysels   PetscBool         iascii;
361ad0c5e61SPieter Ghysels   PetscViewerFormat format;
362ad0c5e61SPieter Ghysels 
363ad0c5e61SPieter Ghysels   PetscFunctionBegin;
364ad0c5e61SPieter Ghysels   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
365ad0c5e61SPieter Ghysels   if (iascii) {
366ad0c5e61SPieter Ghysels     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
367ad0c5e61SPieter Ghysels     if (format == PETSC_VIEWER_ASCII_INFO) {
368860c79edSBarry Smith       ierr = MatView_Info_STRUMPACK(A,viewer);CHKERRQ(ierr);
369ad0c5e61SPieter Ghysels     }
370ad0c5e61SPieter Ghysels   }
371ad0c5e61SPieter Ghysels   PetscFunctionReturn(0);
372ad0c5e61SPieter Ghysels }
3737d6ea485SPieter Ghysels 
374ad0c5e61SPieter Ghysels static PetscErrorCode MatLUFactorNumeric_STRUMPACK(Mat F,Mat A,const MatFactorInfo *info)
3757d6ea485SPieter Ghysels {
3760d08a34dSPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
3777d6ea485SPieter Ghysels   STRUMPACK_RETURN_CODE  sp_err;
3787d6ea485SPieter Ghysels   Mat_SeqAIJ             *A_d,*A_o;
3797d6ea485SPieter Ghysels   Mat_MPIAIJ             *mat;
3807d6ea485SPieter Ghysels   PetscErrorCode         ierr;
3810d08a34dSPieter Ghysels   PetscInt               M=A->rmap->N,m=A->rmap->n;
3827d6ea485SPieter Ghysels   PetscBool              flg;
3837d6ea485SPieter Ghysels 
3847d6ea485SPieter Ghysels   PetscFunctionBegin;
3857d6ea485SPieter Ghysels   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
3867d6ea485SPieter Ghysels   if (flg) { /* A is MATMPIAIJ */
3877d6ea485SPieter Ghysels     mat = (Mat_MPIAIJ*)A->data;
3887d6ea485SPieter Ghysels     A_d = (Mat_SeqAIJ*)(mat->A)->data;
3897d6ea485SPieter Ghysels     A_o = (Mat_SeqAIJ*)(mat->B)->data;
3900d08a34dSPieter 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));
3917d6ea485SPieter Ghysels   } else { /* A is MATSEQAIJ */
3927d6ea485SPieter Ghysels     A_d = (Mat_SeqAIJ*)A->data;
3930d08a34dSPieter Ghysels     PetscStackCall("STRUMPACK_set_csr_matrix",STRUMPACK_set_csr_matrix(*S,&M,A_d->i,A_d->j,A_d->a,0));
3947d6ea485SPieter Ghysels   }
3957d6ea485SPieter Ghysels 
3967d6ea485SPieter Ghysels   /* Reorder and Factor the matrix. */
3977d6ea485SPieter Ghysels   /* TODO figure out how to avoid reorder if the matrix values changed, but the pattern remains the same. */
3980d08a34dSPieter Ghysels   PetscStackCall("STRUMPACK_reorder",sp_err = STRUMPACK_reorder(*S));
3990d08a34dSPieter Ghysels   PetscStackCall("STRUMPACK_factor",sp_err = STRUMPACK_factor(*S));
4000d08a34dSPieter Ghysels   switch (sp_err) {
4010d08a34dSPieter Ghysels   case STRUMPACK_SUCCESS: break;
4020d08a34dSPieter Ghysels   case STRUMPACK_MATRIX_NOT_SET:   { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); break; }
4030d08a34dSPieter Ghysels   case STRUMPACK_REORDERING_ERROR: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); break; }
4040d08a34dSPieter Ghysels   default:                           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: factorization failed");
4057d6ea485SPieter Ghysels   }
406cb250fa3SPieter Ghysels   F->assembled    = PETSC_TRUE;
407cb250fa3SPieter Ghysels   F->preallocated = PETSC_TRUE;
4087d6ea485SPieter Ghysels   PetscFunctionReturn(0);
4097d6ea485SPieter Ghysels }
4107d6ea485SPieter Ghysels 
411ad0c5e61SPieter Ghysels static PetscErrorCode MatLUFactorSymbolic_STRUMPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
4127d6ea485SPieter Ghysels {
4137d6ea485SPieter Ghysels   PetscFunctionBegin;
4147d6ea485SPieter Ghysels   F->ops->lufactornumeric = MatLUFactorNumeric_STRUMPACK;
4157d6ea485SPieter Ghysels   F->ops->solve           = MatSolve_STRUMPACK;
4167d6ea485SPieter Ghysels   F->ops->matsolve        = MatMatSolve_STRUMPACK;
4177d6ea485SPieter Ghysels   PetscFunctionReturn(0);
4187d6ea485SPieter Ghysels }
4197d6ea485SPieter Ghysels 
420ea799195SBarry Smith static PetscErrorCode MatFactorGetSolverType_aij_strumpack(Mat A,MatSolverType *type)
4217d6ea485SPieter Ghysels {
4227d6ea485SPieter Ghysels   PetscFunctionBegin;
4237d6ea485SPieter Ghysels   *type = MATSOLVERSTRUMPACK;
4247d6ea485SPieter Ghysels   PetscFunctionReturn(0);
4257d6ea485SPieter Ghysels }
4267d6ea485SPieter Ghysels 
427575704cbSPieter Ghysels /*MC
428575704cbSPieter Ghysels   MATSOLVERSSTRUMPACK = "strumpack" - A solver package providing a direct sparse solver (PCLU)
429575704cbSPieter Ghysels   and a preconditioner (PCILU) using low-rank compression via the external package STRUMPACK.
430575704cbSPieter Ghysels 
431575704cbSPieter Ghysels   Consult the STRUMPACK-sparse manual for more info.
432575704cbSPieter Ghysels 
433575704cbSPieter Ghysels   Use
434575704cbSPieter Ghysels      ./configure --download-strumpack
435575704cbSPieter Ghysels   to have PETSc installed with STRUMPACK
436575704cbSPieter Ghysels 
437575704cbSPieter Ghysels   Use
4383ca39a21SBarry Smith     -pc_type lu -pc_factor_mat_solver_type strumpack
439575704cbSPieter Ghysels   to use this as an exact (direct) solver, use
4403ca39a21SBarry Smith     -pc_type ilu -pc_factor_mat_solver_type strumpack
441575704cbSPieter Ghysels   to enable low-rank compression (i.e, use as a preconditioner).
442575704cbSPieter Ghysels 
443575704cbSPieter Ghysels   Works with AIJ matrices
444575704cbSPieter Ghysels 
445575704cbSPieter Ghysels   Options Database Keys:
44634f43fa5SPieter Ghysels + -mat_strumpack_verbose
44734f43fa5SPieter Ghysels . -mat_strumpack_hss_rel_tol <1e-2>         - Relative compression tolerance (None)
44834f43fa5SPieter Ghysels . -mat_strumpack_hss_abs_tol <1e-8>         - Absolute compression tolerance (None)
449575704cbSPieter Ghysels . -mat_strumpack_colperm <TRUE>             - Permute matrix to make diagonal nonzeros (None)
45034f43fa5SPieter Ghysels . -mat_strumpack_hss_min_sep_size <256>     - Minimum size of separator for HSS compression (None)
45134f43fa5SPieter Ghysels . -mat_strumpack_max_rank                   - Maximum rank in HSS compression, when using pctype ilu (None)
452a36bf211SPieter Ghysels . -mat_strumpack_leaf_size                  - Size of diagonal blocks in HSS approximation, when using pctype ilu (None)
45334f43fa5SPieter Ghysels . -mat_strumpack_reordering <METIS>         - Sparsity reducing matrix reordering (choose one of) NATURAL METIS PARMETIS SCOTCH PTSCOTCH RCM (None)
454a2b725a8SWilliam Gropp - -mat_strumpack_iterative_solver <DIRECT>  - Select iterative solver from STRUMPACK (choose one of) AUTO DIRECT REFINE PREC_GMRES GMRES PREC_BICGSTAB BICGSTAB (None)
455575704cbSPieter Ghysels 
456575704cbSPieter Ghysels  Level: beginner
457575704cbSPieter Ghysels 
4583ca39a21SBarry Smith .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, PCFactorSetMatSolverType(), MatSolverType
459575704cbSPieter Ghysels M*/
460ad0c5e61SPieter Ghysels static PetscErrorCode MatGetFactor_aij_strumpack(Mat A,MatFactorType ftype,Mat *F)
4617d6ea485SPieter Ghysels {
4627d6ea485SPieter Ghysels   Mat                           B;
4637d6ea485SPieter Ghysels   PetscErrorCode                ierr;
4647d6ea485SPieter Ghysels   PetscInt                      M=A->rmap->N,N=A->cmap->N;
465575704cbSPieter Ghysels   PetscBool                     verb,flg,set;
46634f43fa5SPieter Ghysels   PetscReal                     ctol;
467a36bf211SPieter Ghysels   PetscInt                      hssminsize,max_rank,leaf_size;
46834f43fa5SPieter Ghysels   STRUMPACK_SparseSolver        *S;
46934f43fa5SPieter Ghysels   STRUMPACK_INTERFACE           iface;
47034f43fa5SPieter Ghysels   STRUMPACK_REORDERING_STRATEGY ndcurrent,ndvalue;
47134f43fa5SPieter Ghysels   STRUMPACK_KRYLOV_SOLVER       itcurrent,itsolver;
47234f43fa5SPieter Ghysels   const STRUMPACK_PRECISION     table[2][2][2] =
47334f43fa5SPieter Ghysels     {{{STRUMPACK_FLOATCOMPLEX_64, STRUMPACK_DOUBLECOMPLEX_64},
4747d6ea485SPieter Ghysels       {STRUMPACK_FLOAT_64,        STRUMPACK_DOUBLE_64}},
4757d6ea485SPieter Ghysels      {{STRUMPACK_FLOATCOMPLEX,    STRUMPACK_DOUBLECOMPLEX},
4767d6ea485SPieter Ghysels       {STRUMPACK_FLOAT,           STRUMPACK_DOUBLE}}};
4774ac6704cSBarry Smith   const STRUMPACK_PRECISION     prec = table[(sizeof(PetscInt)==8)?0:1][(PETSC_SCALAR==PETSC_COMPLEX)?0:1][(PETSC_REAL==PETSC_FLOAT)?0:1];
4784ac6704cSBarry Smith   const char *const             STRUMPACKNDTypes[] = {"NATURAL","METIS","PARMETIS","SCOTCH","PTSCOTCH","RCM","STRUMPACKNDTypes","",0};
4794ac6704cSBarry Smith   const char *const             SolverTypes[] = {"AUTO","NONE","REFINE","PREC_GMRES","GMRES","PREC_BICGSTAB","BICGSTAB","SolverTypes","",0};
4807d6ea485SPieter Ghysels 
4817d6ea485SPieter Ghysels   PetscFunctionBegin;
4827d6ea485SPieter Ghysels   /* Create the factorization matrix */
4837d6ea485SPieter Ghysels   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
4847d6ea485SPieter Ghysels   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,M,N);CHKERRQ(ierr);
4857d6ea485SPieter Ghysels   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
4861e1ea65dSPierre Jolivet   ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr);
4877d6ea485SPieter Ghysels   ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr);
48866e17bc3SBarry Smith   B->trivialsymbolic = PETSC_TRUE;
489575704cbSPieter Ghysels   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) {
4907d6ea485SPieter Ghysels     B->ops->lufactorsymbolic  = MatLUFactorSymbolic_STRUMPACK;
491575704cbSPieter Ghysels     B->ops->ilufactorsymbolic = MatLUFactorSymbolic_STRUMPACK;
492575704cbSPieter Ghysels   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
4937d6ea485SPieter Ghysels   B->ops->view        = MatView_STRUMPACK;
4947d6ea485SPieter Ghysels   B->ops->destroy     = MatDestroy_STRUMPACK;
4957d6ea485SPieter Ghysels   B->ops->getdiagonal = MatGetDiagonal_STRUMPACK;
4963ca39a21SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_aij_strumpack);CHKERRQ(ierr);
49734f43fa5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetReordering_C",MatSTRUMPACKSetReordering_STRUMPACK);CHKERRQ(ierr);
49834f43fa5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetColPerm_C",MatSTRUMPACKSetColPerm_STRUMPACK);CHKERRQ(ierr);
49934f43fa5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSRelTol_C",MatSTRUMPACKSetHSSRelTol_STRUMPACK);CHKERRQ(ierr);
50034f43fa5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSAbsTol_C",MatSTRUMPACKSetHSSAbsTol_STRUMPACK);CHKERRQ(ierr);
501a36bf211SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSMaxRank_C",MatSTRUMPACKSetHSSMaxRank_STRUMPACK);CHKERRQ(ierr);
502a36bf211SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSLeafSize_C",MatSTRUMPACKSetHSSLeafSize_STRUMPACK);CHKERRQ(ierr);
503291d0ed5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSMinSepSize_C",MatSTRUMPACKSetHSSMinSepSize_STRUMPACK);CHKERRQ(ierr);
504575704cbSPieter Ghysels   B->factortype = ftype;
5050d08a34dSPieter Ghysels   ierr     = PetscNewLog(B,&S);CHKERRQ(ierr);
5060d08a34dSPieter Ghysels   B->spptr = S;
5070d08a34dSPieter Ghysels 
5081e1ea65dSPierre Jolivet   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr);
5090d08a34dSPieter Ghysels   iface = flg ? STRUMPACK_MT : STRUMPACK_MPI_DIST;
5107d6ea485SPieter Ghysels 
5117d6ea485SPieter Ghysels   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"STRUMPACK Options","Mat");CHKERRQ(ierr);
5127d6ea485SPieter Ghysels 
513575704cbSPieter Ghysels   verb = PetscLogPrintInfo ? PETSC_TRUE : PETSC_FALSE;
5147d6ea485SPieter Ghysels   ierr = PetscOptionsBool("-mat_strumpack_verbose","Print STRUMPACK information","None",verb,&verb,NULL);CHKERRQ(ierr);
5157d6ea485SPieter Ghysels 
51634f43fa5SPieter Ghysels   PetscStackCall("STRUMPACK_init",STRUMPACK_init(S,PetscObjectComm((PetscObject)A),prec,iface,0,NULL,verb));
51755c022e5SPieter Ghysels 
518a36bf211SPieter Ghysels   PetscStackCall("STRUMPACK_HSS_rel_tol",ctol = (PetscReal)STRUMPACK_HSS_rel_tol(*S));
51934f43fa5SPieter Ghysels   ierr = PetscOptionsReal("-mat_strumpack_hss_rel_tol","Relative compression tolerance","None",ctol,&ctol,&set);CHKERRQ(ierr);
520a36bf211SPieter Ghysels   if (set) PetscStackCall("STRUMPACK_set_HSS_rel_tol",STRUMPACK_set_HSS_rel_tol(*S,(double)ctol));
5217d6ea485SPieter Ghysels 
522a36bf211SPieter Ghysels   PetscStackCall("STRUMPACK_HSS_abs_tol",ctol = (PetscReal)STRUMPACK_HSS_abs_tol(*S));
52334f43fa5SPieter Ghysels   ierr = PetscOptionsReal("-mat_strumpack_hss_abs_tol","Absolute compression tolerance","None",ctol,&ctol,&set);CHKERRQ(ierr);
524a36bf211SPieter Ghysels   if (set) PetscStackCall("STRUMPACK_set_HSS_abs_tol",STRUMPACK_set_HSS_abs_tol(*S,(double)ctol));
525575704cbSPieter Ghysels 
526291d0ed5SPieter Ghysels   PetscStackCall("STRUMPACK_mc64job",flg = (STRUMPACK_mc64job(*S) == 0) ? PETSC_FALSE : PETSC_TRUE);
527575704cbSPieter Ghysels   ierr = PetscOptionsBool("-mat_strumpack_colperm","Find a col perm to get nonzero diagonal","None",flg,&flg,&set);CHKERRQ(ierr);
5280d08a34dSPieter Ghysels   if (set) PetscStackCall("STRUMPACK_set_mc64job",STRUMPACK_set_mc64job(*S,flg ? 5 : 0));
529575704cbSPieter Ghysels 
530291d0ed5SPieter Ghysels   PetscStackCall("STRUMPACK_HSS_min_sep_size",hssminsize = (PetscInt)STRUMPACK_HSS_min_sep_size(*S));
531291d0ed5SPieter Ghysels   ierr = PetscOptionsInt("-mat_strumpack_hss_min_sep_size","Minimum size of separator for HSS compression","None",hssminsize,&hssminsize,&set);CHKERRQ(ierr);
532291d0ed5SPieter Ghysels   if (set) PetscStackCall("STRUMPACK_set_HSS_min_sep_size",STRUMPACK_set_HSS_min_sep_size(*S,(int)hssminsize));
533575704cbSPieter Ghysels 
534a36bf211SPieter Ghysels   PetscStackCall("STRUMPACK_HSS_max_rank",max_rank = (PetscInt)STRUMPACK_HSS_max_rank(*S));
53534f43fa5SPieter Ghysels   ierr = PetscOptionsInt("-mat_strumpack_max_rank","Maximum rank in HSS approximation","None",max_rank,&max_rank,&set);CHKERRQ(ierr);
536a36bf211SPieter Ghysels   if (set) PetscStackCall("STRUMPACK_set_HSS_max_rank",STRUMPACK_set_HSS_max_rank(*S,(int)max_rank));
537a36bf211SPieter Ghysels 
538a36bf211SPieter Ghysels   PetscStackCall("STRUMPACK_HSS_leaf_size",leaf_size = (PetscInt)STRUMPACK_HSS_leaf_size(*S));
539a36bf211SPieter Ghysels   ierr = PetscOptionsInt("-mat_strumpack_leaf_size","Size of diagonal blocks in HSS approximation","None",leaf_size,&leaf_size,&set);CHKERRQ(ierr);
540a36bf211SPieter Ghysels   if (set) PetscStackCall("STRUMPACK_set_HSS_leaf_size",STRUMPACK_set_HSS_leaf_size(*S,(int)leaf_size));
54134f43fa5SPieter Ghysels 
54234f43fa5SPieter Ghysels   PetscStackCall("STRUMPACK_reordering_method",ndcurrent = STRUMPACK_reordering_method(*S));
543*2f613bf5SBarry Smith   ierr = PetscOptionsEnum("-mat_strumpack_reordering","Sparsity reducing matrix reordering","None",STRUMPACKNDTypes,(PetscEnum)ndcurrent,(PetscEnum*)&ndvalue,&set);CHKERRQ(ierr);
54434f43fa5SPieter Ghysels   if (set) PetscStackCall("STRUMPACK_set_reordering_method",STRUMPACK_set_reordering_method(*S,ndvalue));
545575704cbSPieter Ghysels 
54688382b05SPieter Ghysels   if (ftype == MAT_FACTOR_ILU) {
54788382b05SPieter Ghysels     /* When enabling HSS compression, the STRUMPACK solver becomes an incomplete              */
54888382b05SPieter Ghysels     /* (or approximate) LU factorization.                                                     */
549291d0ed5SPieter Ghysels     PetscStackCall("STRUMPACK_enable_HSS",STRUMPACK_enable_HSS(*S));
55088382b05SPieter Ghysels   }
551575704cbSPieter Ghysels 
55234f43fa5SPieter Ghysels   /* Disable the outer iterative solver from STRUMPACK.                                       */
55334f43fa5SPieter Ghysels   /* When STRUMPACK is used as a direct solver, it will by default do iterative refinement.   */
55434f43fa5SPieter Ghysels   /* When STRUMPACK is used as an approximate factorization preconditioner (by enabling       */
55534f43fa5SPieter Ghysels   /* low-rank compression), it will use it's own preconditioned GMRES. Here we can disable    */
55634f43fa5SPieter Ghysels   /* the outer iterative solver, as PETSc uses STRUMPACK from within a KSP.                   */
55734f43fa5SPieter Ghysels   PetscStackCall("STRUMPACK_set_Krylov_solver", STRUMPACK_set_Krylov_solver(*S, STRUMPACK_DIRECT));
55834f43fa5SPieter Ghysels 
55934f43fa5SPieter Ghysels   PetscStackCall("STRUMPACK_Krylov_solver",itcurrent = STRUMPACK_Krylov_solver(*S));
560*2f613bf5SBarry Smith   ierr = PetscOptionsEnum("-mat_strumpack_iterative_solver","Select iterative solver from STRUMPACK","None",SolverTypes,(PetscEnum)itcurrent,(PetscEnum*)&itsolver,&set);CHKERRQ(ierr);
56134f43fa5SPieter Ghysels   if (set) PetscStackCall("STRUMPACK_set_Krylov_solver",STRUMPACK_set_Krylov_solver(*S,itsolver));
56234f43fa5SPieter Ghysels 
56334f43fa5SPieter Ghysels   PetscOptionsEnd();
56455c022e5SPieter Ghysels 
5657d6ea485SPieter Ghysels   *F = B;
5667d6ea485SPieter Ghysels   PetscFunctionReturn(0);
5677d6ea485SPieter Ghysels }
5687d6ea485SPieter Ghysels 
5693ca39a21SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_STRUMPACK(void)
5707d6ea485SPieter Ghysels {
5717d6ea485SPieter Ghysels   PetscErrorCode ierr;
5727d6ea485SPieter Ghysels 
5737d6ea485SPieter Ghysels   PetscFunctionBegin;
5743ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr);
5753ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr);
5763ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_ILU,MatGetFactor_aij_strumpack);CHKERRQ(ierr);
5773ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_ILU,MatGetFactor_aij_strumpack);CHKERRQ(ierr);
5787d6ea485SPieter Ghysels   PetscFunctionReturn(0);
5797d6ea485SPieter Ghysels }
580