xref: /petsc/src/mat/impls/aij/mpi/strumpack/strumpack.c (revision 542aee0f1016d40cd405452dff64faed51f7f209)
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 
5291d0ed5SPieter Ghysels #undef __FUNCT__
6291d0ed5SPieter Ghysels #define __FUNCT__ "MatGetDiagonal_STRUMPACK"
7ad0c5e61SPieter Ghysels static PetscErrorCode MatGetDiagonal_STRUMPACK(Mat A,Vec v)
87d6ea485SPieter Ghysels {
97d6ea485SPieter Ghysels   PetscFunctionBegin;
107d6ea485SPieter Ghysels   SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type: STRUMPACK factor");
117d6ea485SPieter Ghysels   PetscFunctionReturn(0);
127d6ea485SPieter Ghysels }
137d6ea485SPieter Ghysels 
14291d0ed5SPieter Ghysels #undef __FUNCT__
15291d0ed5SPieter Ghysels #define __FUNCT__ "MatDestroy_STRUMPACK"
16ad0c5e61SPieter Ghysels static PetscErrorCode MatDestroy_STRUMPACK(Mat A)
177d6ea485SPieter Ghysels {
180d08a34dSPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)A->spptr;
197d6ea485SPieter Ghysels   PetscErrorCode         ierr;
207d6ea485SPieter Ghysels   PetscBool              flg;
217d6ea485SPieter Ghysels 
227d6ea485SPieter Ghysels   PetscFunctionBegin;
237d6ea485SPieter Ghysels   /* Deallocate STRUMPACK storage */
240d08a34dSPieter Ghysels   PetscStackCall("STRUMPACK_destroy",STRUMPACK_destroy(S));
257d6ea485SPieter Ghysels   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
267d6ea485SPieter Ghysels   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr);
277d6ea485SPieter Ghysels   if (flg) {
287d6ea485SPieter Ghysels     ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
297d6ea485SPieter Ghysels   } else {
307d6ea485SPieter Ghysels     ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr);
317d6ea485SPieter Ghysels   }
32575704cbSPieter Ghysels 
33575704cbSPieter Ghysels   /* clear composed functions */
34575704cbSPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr);
3534f43fa5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetReordering_C",NULL);CHKERRQ(ierr);
3634f43fa5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetColPerm_C",NULL);CHKERRQ(ierr);
3734f43fa5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSRelTol_C",NULL);CHKERRQ(ierr);
3834f43fa5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSAbsTol_C",NULL);CHKERRQ(ierr);
3934f43fa5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSMaxRank_C",NULL);CHKERRQ(ierr);
40a36bf211SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSLeafSize_C",NULL);CHKERRQ(ierr);
41291d0ed5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSMinSepSize_C",NULL);CHKERRQ(ierr);
42575704cbSPieter Ghysels 
43575704cbSPieter Ghysels   PetscFunctionReturn(0);
44575704cbSPieter Ghysels }
45575704cbSPieter Ghysels 
4634f43fa5SPieter Ghysels 
47291d0ed5SPieter Ghysels #undef __FUNCT__
4834f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetReordering_STRUMPACK"
4934f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetReordering_STRUMPACK(Mat F,MatSTRUMPACKReordering reordering)
50575704cbSPieter Ghysels {
510d08a34dSPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
52575704cbSPieter Ghysels 
53575704cbSPieter Ghysels   PetscFunctionBegin;
5434f43fa5SPieter Ghysels   PetscStackCall("STRUMPACK_reordering_method",STRUMPACK_set_reordering_method(*S,(STRUMPACK_REORDERING_STRATEGY)reordering));
5534f43fa5SPieter Ghysels   PetscFunctionReturn(0);
5634f43fa5SPieter Ghysels }
5734f43fa5SPieter Ghysels #undef __FUNCT__
5834f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetReordering"
59*542aee0fSPieter Ghysels /*@
6034f43fa5SPieter Ghysels   MatSTRUMPACKSetReordering - Set STRUMPACK fill-reducing reordering
6134f43fa5SPieter Ghysels 
6234f43fa5SPieter Ghysels    Input Parameters:
6334f43fa5SPieter Ghysels +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
6434f43fa5SPieter Ghysels -  reordering - the code to be used to find the fill-reducing reordering
6534f43fa5SPieter Ghysels       Possible values: NATURAL=0 METIS=1 PARMETIS=2 SCOTCH=3 PTSCOTCH=4 RCM=5
6634f43fa5SPieter Ghysels 
6734f43fa5SPieter Ghysels   Options Database:
6834f43fa5SPieter Ghysels .   -mat_strumpack_reordering <METIS>  - Sparsity reducing matrix reordering (choose one of) NATURAL METIS PARMETIS SCOTCH PTSCOTCH RCM (None)
6934f43fa5SPieter Ghysels 
7034f43fa5SPieter Ghysels    Level: beginner
7134f43fa5SPieter Ghysels 
7234f43fa5SPieter Ghysels    References:
7334f43fa5SPieter Ghysels .      STRUMPACK manual
7434f43fa5SPieter Ghysels 
7534f43fa5SPieter Ghysels .seealso: MatGetFactor()
76*542aee0fSPieter Ghysels @*/
7734f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetReordering(Mat F,MatSTRUMPACKReordering reordering)
7834f43fa5SPieter Ghysels {
7934f43fa5SPieter Ghysels   PetscErrorCode ierr;
8034f43fa5SPieter Ghysels 
8134f43fa5SPieter Ghysels   PetscFunctionBegin;
8234f43fa5SPieter Ghysels   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
8334f43fa5SPieter Ghysels   PetscValidLogicalCollectiveEnum(F,reordering,2);
8434f43fa5SPieter Ghysels   ierr = PetscTryMethod(F,"MatSTRUMPACKSetReordering_C",(Mat,MatSTRUMPACKReordering),(F,reordering));CHKERRQ(ierr);
85575704cbSPieter Ghysels   PetscFunctionReturn(0);
86575704cbSPieter Ghysels }
87575704cbSPieter Ghysels 
8834f43fa5SPieter Ghysels 
89291d0ed5SPieter Ghysels #undef __FUNCT__
9034f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetColPerm_STRUMPACK"
9134f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetColPerm_STRUMPACK(Mat F,PetscBool cperm)
9234f43fa5SPieter Ghysels {
9334f43fa5SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
9434f43fa5SPieter Ghysels 
9534f43fa5SPieter Ghysels   PetscFunctionBegin;
9634f43fa5SPieter Ghysels   PetscStackCall("STRUMPACK_set_mc64job",STRUMPACK_set_mc64job(*S,cperm ? 5 : 0));
9734f43fa5SPieter Ghysels   PetscFunctionReturn(0);
9834f43fa5SPieter Ghysels }
9934f43fa5SPieter Ghysels #undef __FUNCT__
10034f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetColPerm"
101575704cbSPieter Ghysels /*@
10234f43fa5SPieter Ghysels   MatSTRUMPACKSetColPerm - Set whether STRUMPACK should try to permute the columns of the matrix in order to get a nonzero diagonal
103575704cbSPieter Ghysels    Logically Collective on Mat
104575704cbSPieter Ghysels 
105575704cbSPieter Ghysels    Input Parameters:
106575704cbSPieter Ghysels +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
10734f43fa5SPieter Ghysels -  cperm - whether or not to permute (internally) the columns of the matrix
108575704cbSPieter Ghysels 
109575704cbSPieter Ghysels   Options Database:
11034f43fa5SPieter Ghysels .   -mat_strumpack_colperm <cperm>
111575704cbSPieter Ghysels 
112575704cbSPieter Ghysels    Level: beginner
113575704cbSPieter Ghysels 
114575704cbSPieter Ghysels    References:
115575704cbSPieter Ghysels .      STRUMPACK manual
116575704cbSPieter Ghysels 
117575704cbSPieter Ghysels .seealso: MatGetFactor()
118575704cbSPieter Ghysels @*/
11934f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetColPerm(Mat F,PetscBool cperm)
120575704cbSPieter Ghysels {
121575704cbSPieter Ghysels   PetscErrorCode ierr;
122575704cbSPieter Ghysels 
123575704cbSPieter Ghysels   PetscFunctionBegin;
124575704cbSPieter Ghysels   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
12534f43fa5SPieter Ghysels   PetscValidLogicalCollectiveBool(F,cperm,2);
12634f43fa5SPieter Ghysels   ierr = PetscTryMethod(F,"MatSTRUMPACKSetColPerm_C",(Mat,PetscBool),(F,cperm));CHKERRQ(ierr);
127575704cbSPieter Ghysels   PetscFunctionReturn(0);
128575704cbSPieter Ghysels }
129575704cbSPieter Ghysels 
130291d0ed5SPieter Ghysels #undef __FUNCT__
13134f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSRelTol_STRUMPACK"
13234f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSRelTol_STRUMPACK(Mat F,PetscReal rtol)
13334f43fa5SPieter Ghysels {
13434f43fa5SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
13534f43fa5SPieter Ghysels 
13634f43fa5SPieter Ghysels   PetscFunctionBegin;
137a36bf211SPieter Ghysels   PetscStackCall("STRUMPACK_set_HSS_rel_tol", STRUMPACK_set_HSS_rel_tol(*S,rtol));
13834f43fa5SPieter Ghysels   PetscFunctionReturn(0);
13934f43fa5SPieter Ghysels }
14034f43fa5SPieter Ghysels #undef __FUNCT__
14134f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSRelTol"
14234f43fa5SPieter Ghysels /*@
14334f43fa5SPieter Ghysels   MatSTRUMPACKSetHSSRelTol - Set STRUMPACK relative tolerance for HSS compression
14434f43fa5SPieter Ghysels    Logically Collective on Mat
14534f43fa5SPieter Ghysels 
14634f43fa5SPieter Ghysels    Input Parameters:
14734f43fa5SPieter Ghysels +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
14834f43fa5SPieter Ghysels -  rtol - relative compression tolerance
14934f43fa5SPieter Ghysels 
15034f43fa5SPieter Ghysels   Options Database:
15134f43fa5SPieter Ghysels .   -mat_strumpack_hss_rel_tol <1e-2>         - Relative compression tolerance (None)
15234f43fa5SPieter Ghysels 
15334f43fa5SPieter Ghysels    Level: beginner
15434f43fa5SPieter Ghysels 
15534f43fa5SPieter Ghysels    References:
15634f43fa5SPieter Ghysels .      STRUMPACK manual
15734f43fa5SPieter Ghysels 
15834f43fa5SPieter Ghysels .seealso: MatGetFactor()
15934f43fa5SPieter Ghysels @*/
16034f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSRelTol(Mat F,PetscReal rtol)
16134f43fa5SPieter Ghysels {
16234f43fa5SPieter Ghysels   PetscErrorCode ierr;
16334f43fa5SPieter Ghysels 
16434f43fa5SPieter Ghysels   PetscFunctionBegin;
16534f43fa5SPieter Ghysels   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
16634f43fa5SPieter Ghysels   PetscValidLogicalCollectiveReal(F,rtol,2);
16734f43fa5SPieter Ghysels   ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSRelTol_C",(Mat,PetscReal),(F,rtol));CHKERRQ(ierr);
16834f43fa5SPieter Ghysels   PetscFunctionReturn(0);
16934f43fa5SPieter Ghysels }
17034f43fa5SPieter Ghysels 
17134f43fa5SPieter Ghysels 
17234f43fa5SPieter Ghysels #undef __FUNCT__
17334f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSAbsTol_STRUMPACK"
17434f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSAbsTol_STRUMPACK(Mat F,PetscReal atol)
17534f43fa5SPieter Ghysels {
17634f43fa5SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
17734f43fa5SPieter Ghysels 
17834f43fa5SPieter Ghysels   PetscFunctionBegin;
179a36bf211SPieter Ghysels   PetscStackCall("STRUMPACK_set_HSS_abs_tol", STRUMPACK_set_HSS_abs_tol(*S,atol));
18034f43fa5SPieter Ghysels   PetscFunctionReturn(0);
18134f43fa5SPieter Ghysels }
18234f43fa5SPieter Ghysels #undef __FUNCT__
18334f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSAbsTol"
18434f43fa5SPieter Ghysels /*@
18534f43fa5SPieter Ghysels   MatSTRUMPACKSetHSSAbsTol - Set STRUMPACK absolute tolerance for HSS compression
18634f43fa5SPieter Ghysels    Logically Collective on Mat
18734f43fa5SPieter Ghysels 
18834f43fa5SPieter Ghysels    Input Parameters:
18934f43fa5SPieter Ghysels +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
19034f43fa5SPieter Ghysels -  atol - absolute compression tolerance
19134f43fa5SPieter Ghysels 
19234f43fa5SPieter Ghysels   Options Database:
19334f43fa5SPieter Ghysels .   -mat_strumpack_hss_abs_tol <1e-8>         - Absolute compression tolerance (None)
19434f43fa5SPieter Ghysels 
19534f43fa5SPieter Ghysels    Level: beginner
19634f43fa5SPieter Ghysels 
19734f43fa5SPieter Ghysels    References:
19834f43fa5SPieter Ghysels .      STRUMPACK manual
19934f43fa5SPieter Ghysels 
20034f43fa5SPieter Ghysels .seealso: MatGetFactor()
20134f43fa5SPieter Ghysels @*/
20234f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSAbsTol(Mat F,PetscReal atol)
20334f43fa5SPieter Ghysels {
20434f43fa5SPieter Ghysels   PetscErrorCode ierr;
20534f43fa5SPieter Ghysels 
20634f43fa5SPieter Ghysels   PetscFunctionBegin;
20734f43fa5SPieter Ghysels   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
20834f43fa5SPieter Ghysels   PetscValidLogicalCollectiveReal(F,atol,2);
20934f43fa5SPieter Ghysels   ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSAbsTol_C",(Mat,PetscReal),(F,atol));CHKERRQ(ierr);
21034f43fa5SPieter Ghysels   PetscFunctionReturn(0);
21134f43fa5SPieter Ghysels }
21234f43fa5SPieter Ghysels 
21334f43fa5SPieter Ghysels 
21434f43fa5SPieter Ghysels #undef __FUNCT__
21534f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSMaxRank_STRUMPACK"
21634f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSMaxRank_STRUMPACK(Mat F,PetscInt hssmaxrank)
21734f43fa5SPieter Ghysels {
21834f43fa5SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
21934f43fa5SPieter Ghysels 
22034f43fa5SPieter Ghysels   PetscFunctionBegin;
221a36bf211SPieter Ghysels   PetscStackCall("STRUMPACK_set_HSS_max_rank", STRUMPACK_set_HSS_max_rank(*S,hssmaxrank));
22234f43fa5SPieter Ghysels   PetscFunctionReturn(0);
22334f43fa5SPieter Ghysels }
22434f43fa5SPieter Ghysels #undef __FUNCT__
22534f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSMaxRank"
22634f43fa5SPieter Ghysels /*@
22734f43fa5SPieter Ghysels   MatSTRUMPACKSetHSSMaxRank - Set STRUMPACK maximum HSS rank
22834f43fa5SPieter Ghysels    Logically Collective on Mat
22934f43fa5SPieter Ghysels 
23034f43fa5SPieter Ghysels    Input Parameters:
23134f43fa5SPieter Ghysels +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
23234f43fa5SPieter Ghysels -  hssmaxrank - maximum rank used in low-rank approximation
23334f43fa5SPieter Ghysels 
23434f43fa5SPieter Ghysels   Options Database:
23534f43fa5SPieter Ghysels .   -mat_strumpack_max_rank    - Maximum rank in HSS compression, when using pctype ilu (None)
23634f43fa5SPieter Ghysels 
23734f43fa5SPieter Ghysels    Level: beginner
23834f43fa5SPieter Ghysels 
23934f43fa5SPieter Ghysels    References:
24034f43fa5SPieter Ghysels .      STRUMPACK manual
24134f43fa5SPieter Ghysels 
24234f43fa5SPieter Ghysels .seealso: MatGetFactor()
24334f43fa5SPieter Ghysels @*/
24434f43fa5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSMaxRank(Mat F,PetscInt hssmaxrank)
24534f43fa5SPieter Ghysels {
24634f43fa5SPieter Ghysels   PetscErrorCode ierr;
24734f43fa5SPieter Ghysels 
24834f43fa5SPieter Ghysels   PetscFunctionBegin;
24934f43fa5SPieter Ghysels   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
25034f43fa5SPieter Ghysels   PetscValidLogicalCollectiveInt(F,hssmaxrank,2);
25134f43fa5SPieter Ghysels   ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSMaxRank_C",(Mat,PetscInt),(F,hssmaxrank));CHKERRQ(ierr);
25234f43fa5SPieter Ghysels   PetscFunctionReturn(0);
25334f43fa5SPieter Ghysels }
25434f43fa5SPieter Ghysels 
25534f43fa5SPieter Ghysels 
25634f43fa5SPieter Ghysels #undef __FUNCT__
257a36bf211SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSLeafSize_STRUMPACK"
258a36bf211SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSLeafSize_STRUMPACK(Mat F,PetscInt leaf_size)
259a36bf211SPieter Ghysels {
260a36bf211SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
261a36bf211SPieter Ghysels 
262a36bf211SPieter Ghysels   PetscFunctionBegin;
263a36bf211SPieter Ghysels   PetscStackCall("STRUMPACK_set_HSS_leaf_size", STRUMPACK_set_HSS_leaf_size(*S,leaf_size));
264a36bf211SPieter Ghysels   PetscFunctionReturn(0);
265a36bf211SPieter Ghysels }
266a36bf211SPieter Ghysels #undef __FUNCT__
267a36bf211SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSLeafSize"
268a36bf211SPieter Ghysels /*@
269a36bf211SPieter Ghysels   MatSTRUMPACKSetHSSLeafSize - Set STRUMPACK HSS leaf size
270a36bf211SPieter Ghysels    Logically Collective on Mat
271a36bf211SPieter Ghysels 
272a36bf211SPieter Ghysels    Input Parameters:
273a36bf211SPieter Ghysels +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
274a36bf211SPieter Ghysels -  leaf_size - Size of diagonal blocks in HSS approximation
275a36bf211SPieter Ghysels 
276a36bf211SPieter Ghysels   Options Database:
277a36bf211SPieter Ghysels .   -mat_strumpack_leaf_size    - Size of diagonal blocks in HSS approximation, when using pctype ilu (None)
278a36bf211SPieter Ghysels 
279a36bf211SPieter Ghysels    Level: beginner
280a36bf211SPieter Ghysels 
281a36bf211SPieter Ghysels    References:
282a36bf211SPieter Ghysels .      STRUMPACK manual
283a36bf211SPieter Ghysels 
284a36bf211SPieter Ghysels .seealso: MatGetFactor()
285a36bf211SPieter Ghysels @*/
286a36bf211SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSLeafSize(Mat F,PetscInt leaf_size)
287a36bf211SPieter Ghysels {
288a36bf211SPieter Ghysels   PetscErrorCode ierr;
289a36bf211SPieter Ghysels 
290a36bf211SPieter Ghysels   PetscFunctionBegin;
291a36bf211SPieter Ghysels   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
292a36bf211SPieter Ghysels   PetscValidLogicalCollectiveInt(F,leaf_size,2);
293a36bf211SPieter Ghysels   ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSLeafSize_C",(Mat,PetscInt),(F,leaf_size));CHKERRQ(ierr);
294a36bf211SPieter Ghysels   PetscFunctionReturn(0);
295a36bf211SPieter Ghysels }
296a36bf211SPieter Ghysels 
297a36bf211SPieter Ghysels 
298575704cbSPieter Ghysels 
299291d0ed5SPieter Ghysels #undef __FUNCT__
30034f43fa5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSMinSepSize_STRUMPACK"
30134f43fa5SPieter Ghysels static PetscErrorCode MatSTRUMPACKSetHSSMinSepSize_STRUMPACK(Mat F,PetscInt hssminsize)
30234f43fa5SPieter Ghysels {
30334f43fa5SPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
30434f43fa5SPieter Ghysels 
30534f43fa5SPieter Ghysels   PetscFunctionBegin;
30634f43fa5SPieter Ghysels   PetscStackCall("STRUMPACK_set_HSS_min_sep_size", STRUMPACK_set_HSS_min_sep_size(*S,hssminsize));
30734f43fa5SPieter Ghysels   PetscFunctionReturn(0);
30834f43fa5SPieter Ghysels }
30934f43fa5SPieter Ghysels #undef __FUNCT__
310291d0ed5SPieter Ghysels #define __FUNCT__ "MatSTRUMPACKSetHSSMinSepSize"
311291d0ed5SPieter Ghysels /*@
312291d0ed5SPieter Ghysels   MatSTRUMPACKSetHSSMinSepSize - Set STRUMPACK minimum separator size for low-rank approximation
313291d0ed5SPieter Ghysels    Logically Collective on Mat
314291d0ed5SPieter Ghysels 
315291d0ed5SPieter Ghysels    Input Parameters:
316291d0ed5SPieter Ghysels +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
317291d0ed5SPieter Ghysels -  hssminsize - minimum dense matrix size for low-rank approximation
318291d0ed5SPieter Ghysels 
319291d0ed5SPieter Ghysels   Options Database:
320291d0ed5SPieter Ghysels .   -mat_strumpack_hss_min_sep_size <hssminsize>
321291d0ed5SPieter Ghysels 
322291d0ed5SPieter Ghysels    Level: beginner
323291d0ed5SPieter Ghysels 
324291d0ed5SPieter Ghysels    References:
325291d0ed5SPieter Ghysels .      STRUMPACK manual
326291d0ed5SPieter Ghysels 
327291d0ed5SPieter Ghysels .seealso: MatGetFactor()
328291d0ed5SPieter Ghysels @*/
329291d0ed5SPieter Ghysels PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat F,PetscInt hssminsize)
330291d0ed5SPieter Ghysels {
331291d0ed5SPieter Ghysels   PetscErrorCode ierr;
332291d0ed5SPieter Ghysels 
333291d0ed5SPieter Ghysels   PetscFunctionBegin;
334291d0ed5SPieter Ghysels   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
335291d0ed5SPieter Ghysels   PetscValidLogicalCollectiveInt(F,hssminsize,2);
336291d0ed5SPieter Ghysels   ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSMinSepSize_C",(Mat,PetscInt),(F,hssminsize));CHKERRQ(ierr);
337291d0ed5SPieter Ghysels   PetscFunctionReturn(0);
338291d0ed5SPieter Ghysels }
339291d0ed5SPieter Ghysels 
340291d0ed5SPieter Ghysels 
341291d0ed5SPieter Ghysels #undef __FUNCT__
342291d0ed5SPieter Ghysels #define __FUNCT__ "MatSolve_STRUMPACK"
343ad0c5e61SPieter Ghysels static PetscErrorCode MatSolve_STRUMPACK(Mat A,Vec b_mpi,Vec x)
3447d6ea485SPieter Ghysels {
3450d08a34dSPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)A->spptr;
3467d6ea485SPieter Ghysels   STRUMPACK_RETURN_CODE  sp_err;
3477d6ea485SPieter Ghysels   PetscErrorCode         ierr;
3487d6ea485SPieter Ghysels   const PetscScalar      *bptr;
3497d6ea485SPieter Ghysels   PetscScalar            *xptr;
3507d6ea485SPieter Ghysels 
3517d6ea485SPieter Ghysels   PetscFunctionBegin;
3527d6ea485SPieter Ghysels   ierr = VecGetArray(x,&xptr);CHKERRQ(ierr);
3537d6ea485SPieter Ghysels   ierr = VecGetArrayRead(b_mpi,&bptr);CHKERRQ(ierr);
3540d08a34dSPieter Ghysels 
3550d08a34dSPieter Ghysels   PetscStackCall("STRUMPACK_solve",sp_err = STRUMPACK_solve(*S,(PetscScalar*)bptr,xptr,0));
3560d08a34dSPieter Ghysels   switch (sp_err) {
3570d08a34dSPieter Ghysels   case STRUMPACK_SUCCESS: break;
3580d08a34dSPieter Ghysels   case STRUMPACK_MATRIX_NOT_SET:   { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); break; }
3590d08a34dSPieter Ghysels   case STRUMPACK_REORDERING_ERROR: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); break; }
3600d08a34dSPieter Ghysels   default:                           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: solve failed");
3617d6ea485SPieter Ghysels   }
3627d6ea485SPieter Ghysels   ierr = VecRestoreArray(x,&xptr);CHKERRQ(ierr);
3637d6ea485SPieter Ghysels   ierr = VecRestoreArrayRead(b_mpi,&bptr);CHKERRQ(ierr);
3647d6ea485SPieter Ghysels   PetscFunctionReturn(0);
3657d6ea485SPieter Ghysels }
3667d6ea485SPieter Ghysels 
367291d0ed5SPieter Ghysels #undef __FUNCT__
368291d0ed5SPieter Ghysels #define __FUNCT__ "MatMatSolve_STRUMPACK"
369ad0c5e61SPieter Ghysels static PetscErrorCode MatMatSolve_STRUMPACK(Mat A,Mat B_mpi,Mat X)
3707d6ea485SPieter Ghysels {
3717d6ea485SPieter Ghysels   PetscErrorCode   ierr;
3727d6ea485SPieter Ghysels   PetscBool        flg;
3737d6ea485SPieter Ghysels 
3747d6ea485SPieter Ghysels   PetscFunctionBegin;
3757d6ea485SPieter Ghysels   ierr = PetscObjectTypeCompareAny((PetscObject)B_mpi,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
3767d6ea485SPieter Ghysels   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
3777d6ea485SPieter Ghysels   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
3787d6ea485SPieter Ghysels   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
3797d6ea485SPieter Ghysels   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_STRUMPACK() is not implemented yet");
3807d6ea485SPieter Ghysels   PetscFunctionReturn(0);
3817d6ea485SPieter Ghysels }
3827d6ea485SPieter Ghysels 
383291d0ed5SPieter Ghysels #undef __FUNCT__
384291d0ed5SPieter Ghysels #define __FUNCT__ "MatFactorInfo_STRUMPACK"
385ad0c5e61SPieter Ghysels static PetscErrorCode MatFactorInfo_STRUMPACK(Mat A,PetscViewer viewer)
386ad0c5e61SPieter Ghysels {
387ad0c5e61SPieter Ghysels   PetscErrorCode  ierr;
388ad0c5e61SPieter Ghysels 
389ad0c5e61SPieter Ghysels   PetscFunctionBegin;
390ad0c5e61SPieter Ghysels   /* check if matrix is strumpack type */
391ad0c5e61SPieter Ghysels   if (A->ops->solve != MatSolve_STRUMPACK) PetscFunctionReturn(0);
392ad0c5e61SPieter Ghysels   ierr = PetscViewerASCIIPrintf(viewer,"STRUMPACK sparse solver!\n");CHKERRQ(ierr);
393ad0c5e61SPieter Ghysels   PetscFunctionReturn(0);
394ad0c5e61SPieter Ghysels }
395ad0c5e61SPieter Ghysels 
396291d0ed5SPieter Ghysels #undef __FUNCT__
397291d0ed5SPieter Ghysels #define __FUNCT__ "MatView_STRUMPACK"
398ad0c5e61SPieter Ghysels static PetscErrorCode MatView_STRUMPACK(Mat A,PetscViewer viewer)
399ad0c5e61SPieter Ghysels {
400ad0c5e61SPieter Ghysels   PetscErrorCode    ierr;
401ad0c5e61SPieter Ghysels   PetscBool         iascii;
402ad0c5e61SPieter Ghysels   PetscViewerFormat format;
403ad0c5e61SPieter Ghysels 
404ad0c5e61SPieter Ghysels   PetscFunctionBegin;
405ad0c5e61SPieter Ghysels   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
406ad0c5e61SPieter Ghysels   if (iascii) {
407ad0c5e61SPieter Ghysels     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
408ad0c5e61SPieter Ghysels     if (format == PETSC_VIEWER_ASCII_INFO) {
409ad0c5e61SPieter Ghysels       ierr = MatFactorInfo_STRUMPACK(A,viewer);CHKERRQ(ierr);
410ad0c5e61SPieter Ghysels     }
411ad0c5e61SPieter Ghysels   }
412ad0c5e61SPieter Ghysels   PetscFunctionReturn(0);
413ad0c5e61SPieter Ghysels }
4147d6ea485SPieter Ghysels 
415291d0ed5SPieter Ghysels #undef __FUNCT__
416291d0ed5SPieter Ghysels #define __FUNCT__ "MatLUFactorNumeric_STRUMPACK"
417ad0c5e61SPieter Ghysels static PetscErrorCode MatLUFactorNumeric_STRUMPACK(Mat F,Mat A,const MatFactorInfo *info)
4187d6ea485SPieter Ghysels {
4190d08a34dSPieter Ghysels   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
4207d6ea485SPieter Ghysels   STRUMPACK_RETURN_CODE  sp_err;
4217d6ea485SPieter Ghysels   Mat_SeqAIJ             *A_d,*A_o;
4227d6ea485SPieter Ghysels   Mat_MPIAIJ             *mat;
4237d6ea485SPieter Ghysels   PetscErrorCode         ierr;
4240d08a34dSPieter Ghysels   PetscInt               M=A->rmap->N,m=A->rmap->n;
4257d6ea485SPieter Ghysels   PetscBool              flg;
4267d6ea485SPieter Ghysels 
4277d6ea485SPieter Ghysels   PetscFunctionBegin;
4287d6ea485SPieter Ghysels   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
4297d6ea485SPieter Ghysels   if (flg) { /* A is MATMPIAIJ */
4307d6ea485SPieter Ghysels     mat = (Mat_MPIAIJ*)A->data;
4317d6ea485SPieter Ghysels     A_d = (Mat_SeqAIJ*)(mat->A)->data;
4327d6ea485SPieter Ghysels     A_o = (Mat_SeqAIJ*)(mat->B)->data;
4330d08a34dSPieter 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));
4347d6ea485SPieter Ghysels   } else { /* A is MATSEQAIJ */
4357d6ea485SPieter Ghysels     A_d = (Mat_SeqAIJ*)A->data;
4360d08a34dSPieter Ghysels     PetscStackCall("STRUMPACK_set_csr_matrix",STRUMPACK_set_csr_matrix(*S,&M,A_d->i,A_d->j,A_d->a,0));
4377d6ea485SPieter Ghysels   }
4387d6ea485SPieter Ghysels 
4397d6ea485SPieter Ghysels   /* Reorder and Factor the matrix. */
4407d6ea485SPieter Ghysels   /* TODO figure out how to avoid reorder if the matrix values changed, but the pattern remains the same. */
4410d08a34dSPieter Ghysels   PetscStackCall("STRUMPACK_reorder",sp_err = STRUMPACK_reorder(*S));
4420d08a34dSPieter Ghysels   PetscStackCall("STRUMPACK_factor",sp_err = STRUMPACK_factor(*S));
4430d08a34dSPieter Ghysels   switch (sp_err) {
4440d08a34dSPieter Ghysels   case STRUMPACK_SUCCESS: break;
4450d08a34dSPieter Ghysels   case STRUMPACK_MATRIX_NOT_SET:   { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); break; }
4460d08a34dSPieter Ghysels   case STRUMPACK_REORDERING_ERROR: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); break; }
4470d08a34dSPieter Ghysels   default:                           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: factorization failed");
4487d6ea485SPieter Ghysels   }
4497d6ea485SPieter Ghysels   PetscFunctionReturn(0);
4507d6ea485SPieter Ghysels }
4517d6ea485SPieter Ghysels 
452291d0ed5SPieter Ghysels #undef __FUNCT__
453291d0ed5SPieter Ghysels #define __FUNCT__ "MatLUFactorSymbolic_STRUMPACK"
454ad0c5e61SPieter Ghysels static PetscErrorCode MatLUFactorSymbolic_STRUMPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
4557d6ea485SPieter Ghysels {
4567d6ea485SPieter Ghysels   PetscFunctionBegin;
4577d6ea485SPieter Ghysels   F->ops->lufactornumeric = MatLUFactorNumeric_STRUMPACK;
4587d6ea485SPieter Ghysels   F->ops->solve           = MatSolve_STRUMPACK;
4597d6ea485SPieter Ghysels   F->ops->matsolve        = MatMatSolve_STRUMPACK;
4607d6ea485SPieter Ghysels   PetscFunctionReturn(0);
4617d6ea485SPieter Ghysels }
4627d6ea485SPieter Ghysels 
463291d0ed5SPieter Ghysels #undef __FUNCT__
464291d0ed5SPieter Ghysels #define __FUNCT__ "MatFactorGetSolverPackage_aij_strumpack"
465ad0c5e61SPieter Ghysels static PetscErrorCode MatFactorGetSolverPackage_aij_strumpack(Mat A,const MatSolverPackage *type)
4667d6ea485SPieter Ghysels {
4677d6ea485SPieter Ghysels   PetscFunctionBegin;
4687d6ea485SPieter Ghysels   *type = MATSOLVERSTRUMPACK;
4697d6ea485SPieter Ghysels   PetscFunctionReturn(0);
4707d6ea485SPieter Ghysels }
4717d6ea485SPieter Ghysels 
472575704cbSPieter Ghysels /*MC
473575704cbSPieter Ghysels   MATSOLVERSSTRUMPACK = "strumpack" - A solver package providing a direct sparse solver (PCLU)
474575704cbSPieter Ghysels   and a preconditioner (PCILU) using low-rank compression via the external package STRUMPACK.
475575704cbSPieter Ghysels 
476575704cbSPieter Ghysels   Consult the STRUMPACK-sparse manual for more info.
477575704cbSPieter Ghysels 
478575704cbSPieter Ghysels   Use
479575704cbSPieter Ghysels      ./configure --download-strumpack
480575704cbSPieter Ghysels   to have PETSc installed with STRUMPACK
481575704cbSPieter Ghysels 
482575704cbSPieter Ghysels   Use
483575704cbSPieter Ghysels     -pc_type lu -pc_factor_mat_solver_package strumpack
484575704cbSPieter Ghysels   to use this as an exact (direct) solver, use
485575704cbSPieter Ghysels     -pc_type ilu -pc_factor_mat_solver_package strumpack
486575704cbSPieter Ghysels   to enable low-rank compression (i.e, use as a preconditioner).
487575704cbSPieter Ghysels 
488575704cbSPieter Ghysels   Works with AIJ matrices
489575704cbSPieter Ghysels 
490575704cbSPieter Ghysels   Options Database Keys:
49134f43fa5SPieter Ghysels + -mat_strumpack_verbose
49234f43fa5SPieter Ghysels . -mat_strumpack_hss_rel_tol <1e-2>         - Relative compression tolerance (None)
49334f43fa5SPieter Ghysels . -mat_strumpack_hss_abs_tol <1e-8>         - Absolute compression tolerance (None)
494575704cbSPieter Ghysels . -mat_strumpack_colperm <TRUE>             - Permute matrix to make diagonal nonzeros (None)
49534f43fa5SPieter Ghysels . -mat_strumpack_hss_min_sep_size <256>     - Minimum size of separator for HSS compression (None)
49634f43fa5SPieter Ghysels . -mat_strumpack_max_rank                   - Maximum rank in HSS compression, when using pctype ilu (None)
497a36bf211SPieter Ghysels . -mat_strumpack_leaf_size                  - Size of diagonal blocks in HSS approximation, when using pctype ilu (None)
49834f43fa5SPieter Ghysels . -mat_strumpack_reordering <METIS>         - Sparsity reducing matrix reordering (choose one of) NATURAL METIS PARMETIS SCOTCH PTSCOTCH RCM (None)
49934f43fa5SPieter Ghysels . -mat_strumpack_iterative_solver <DIRECT>  - Select iterative solver from STRUMPACK (choose one of) AUTO DIRECT REFINE PREC_GMRES GMRES PREC_BICGSTAB BICGSTAB (None)
500575704cbSPieter Ghysels 
501575704cbSPieter Ghysels  Level: beginner
502575704cbSPieter Ghysels 
503575704cbSPieter Ghysels .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, PCFactorSetMatSolverPackage(), MatSolverPackage
504575704cbSPieter Ghysels M*/
505291d0ed5SPieter Ghysels #undef __FUNCT__
506291d0ed5SPieter Ghysels #define __FUNCT__ "MatGetFactor_aij_strumpack"
507ad0c5e61SPieter Ghysels static PetscErrorCode MatGetFactor_aij_strumpack(Mat A,MatFactorType ftype,Mat *F)
5087d6ea485SPieter Ghysels {
5097d6ea485SPieter Ghysels   Mat                           B;
5107d6ea485SPieter Ghysels   PetscErrorCode                ierr;
5117d6ea485SPieter Ghysels   PetscInt                      M=A->rmap->N,N=A->cmap->N;
512575704cbSPieter Ghysels   PetscBool                     verb,flg,set;
51334f43fa5SPieter Ghysels   PetscReal                     ctol;
514a36bf211SPieter Ghysels   PetscInt                      hssminsize,max_rank,leaf_size;
51534f43fa5SPieter Ghysels   STRUMPACK_SparseSolver        *S;
51634f43fa5SPieter Ghysels   STRUMPACK_INTERFACE           iface;
51734f43fa5SPieter Ghysels   STRUMPACK_REORDERING_STRATEGY ndcurrent,ndvalue;
51834f43fa5SPieter Ghysels   STRUMPACK_KRYLOV_SOLVER       itcurrent,itsolver;
51934f43fa5SPieter Ghysels     const STRUMPACK_PRECISION table[2][2][2] =
52034f43fa5SPieter Ghysels     {{{STRUMPACK_FLOATCOMPLEX_64, STRUMPACK_DOUBLECOMPLEX_64},
5217d6ea485SPieter Ghysels       {STRUMPACK_FLOAT_64,        STRUMPACK_DOUBLE_64}},
5227d6ea485SPieter Ghysels      {{STRUMPACK_FLOATCOMPLEX,    STRUMPACK_DOUBLECOMPLEX},
5237d6ea485SPieter Ghysels       {STRUMPACK_FLOAT,           STRUMPACK_DOUBLE}}};
52434f43fa5SPieter Ghysels   const STRUMPACK_PRECISION prec =
52534f43fa5SPieter Ghysels     table[(sizeof(PetscInt)==8)?0:1]
52634f43fa5SPieter Ghysels     [(PETSC_SCALAR==PETSC_COMPLEX)?0:1]
52734f43fa5SPieter Ghysels     [(PETSC_REAL==PETSC_FLOAT)?0:1];
5287d6ea485SPieter Ghysels 
5297d6ea485SPieter Ghysels   PetscFunctionBegin;
5307d6ea485SPieter Ghysels   /* Create the factorization matrix */
5317d6ea485SPieter Ghysels   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
5327d6ea485SPieter Ghysels   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,M,N);CHKERRQ(ierr);
5337d6ea485SPieter Ghysels   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
5347d6ea485SPieter Ghysels   ierr = MatSeqAIJSetPreallocation(B,0,NULL);
5357d6ea485SPieter Ghysels   ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr);
536575704cbSPieter Ghysels   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) {
5377d6ea485SPieter Ghysels     B->ops->lufactorsymbolic  = MatLUFactorSymbolic_STRUMPACK;
538575704cbSPieter Ghysels     B->ops->ilufactorsymbolic = MatLUFactorSymbolic_STRUMPACK;
539575704cbSPieter Ghysels   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
5407d6ea485SPieter Ghysels   B->ops->view        = MatView_STRUMPACK;
5417d6ea485SPieter Ghysels   B->ops->destroy     = MatDestroy_STRUMPACK;
5427d6ea485SPieter Ghysels   B->ops->getdiagonal = MatGetDiagonal_STRUMPACK;
5437d6ea485SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_aij_strumpack);CHKERRQ(ierr);
54434f43fa5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetReordering_C",MatSTRUMPACKSetReordering_STRUMPACK);CHKERRQ(ierr);
54534f43fa5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetColPerm_C",MatSTRUMPACKSetColPerm_STRUMPACK);CHKERRQ(ierr);
54634f43fa5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSRelTol_C",MatSTRUMPACKSetHSSRelTol_STRUMPACK);CHKERRQ(ierr);
54734f43fa5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSAbsTol_C",MatSTRUMPACKSetHSSAbsTol_STRUMPACK);CHKERRQ(ierr);
548a36bf211SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSMaxRank_C",MatSTRUMPACKSetHSSMaxRank_STRUMPACK);CHKERRQ(ierr);
549a36bf211SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSLeafSize_C",MatSTRUMPACKSetHSSLeafSize_STRUMPACK);CHKERRQ(ierr);
550291d0ed5SPieter Ghysels   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSMinSepSize_C",MatSTRUMPACKSetHSSMinSepSize_STRUMPACK);CHKERRQ(ierr);
551575704cbSPieter Ghysels   B->factortype = ftype;
5520d08a34dSPieter Ghysels   ierr     = PetscNewLog(B,&S);CHKERRQ(ierr);
5530d08a34dSPieter Ghysels   B->spptr = S;
5540d08a34dSPieter Ghysels 
5550d08a34dSPieter Ghysels   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);
5560d08a34dSPieter Ghysels   iface = flg ? STRUMPACK_MT : STRUMPACK_MPI_DIST;
5577d6ea485SPieter Ghysels 
5587d6ea485SPieter Ghysels   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"STRUMPACK Options","Mat");CHKERRQ(ierr);
5597d6ea485SPieter Ghysels 
560575704cbSPieter Ghysels   verb = PetscLogPrintInfo ? PETSC_TRUE : PETSC_FALSE;
5617d6ea485SPieter Ghysels   ierr = PetscOptionsBool("-mat_strumpack_verbose","Print STRUMPACK information","None",verb,&verb,NULL);CHKERRQ(ierr);
5627d6ea485SPieter Ghysels 
56334f43fa5SPieter Ghysels   PetscStackCall("STRUMPACK_init",STRUMPACK_init(S,PetscObjectComm((PetscObject)A),prec,iface,0,NULL,verb));
56455c022e5SPieter Ghysels 
565a36bf211SPieter Ghysels   PetscStackCall("STRUMPACK_HSS_rel_tol",ctol = (PetscReal)STRUMPACK_HSS_rel_tol(*S));
56634f43fa5SPieter Ghysels   ierr = PetscOptionsReal("-mat_strumpack_hss_rel_tol","Relative compression tolerance","None",ctol,&ctol,&set);CHKERRQ(ierr);
567a36bf211SPieter Ghysels   if (set) PetscStackCall("STRUMPACK_set_HSS_rel_tol",STRUMPACK_set_HSS_rel_tol(*S,(double)ctol));
5687d6ea485SPieter Ghysels 
569a36bf211SPieter Ghysels   PetscStackCall("STRUMPACK_HSS_abs_tol",ctol = (PetscReal)STRUMPACK_HSS_abs_tol(*S));
57034f43fa5SPieter Ghysels   ierr = PetscOptionsReal("-mat_strumpack_hss_abs_tol","Absolute compression tolerance","None",ctol,&ctol,&set);CHKERRQ(ierr);
571a36bf211SPieter Ghysels   if (set) PetscStackCall("STRUMPACK_set_HSS_abs_tol",STRUMPACK_set_HSS_abs_tol(*S,(double)ctol));
572575704cbSPieter Ghysels 
573291d0ed5SPieter Ghysels   PetscStackCall("STRUMPACK_mc64job",flg = (STRUMPACK_mc64job(*S) == 0) ? PETSC_FALSE : PETSC_TRUE);
574575704cbSPieter Ghysels   ierr = PetscOptionsBool("-mat_strumpack_colperm","Find a col perm to get nonzero diagonal","None",flg,&flg,&set);CHKERRQ(ierr);
5750d08a34dSPieter Ghysels   if (set) PetscStackCall("STRUMPACK_set_mc64job",STRUMPACK_set_mc64job(*S,flg ? 5 : 0));
576575704cbSPieter Ghysels 
577291d0ed5SPieter Ghysels   PetscStackCall("STRUMPACK_HSS_min_sep_size",hssminsize = (PetscInt)STRUMPACK_HSS_min_sep_size(*S));
578291d0ed5SPieter Ghysels   ierr = PetscOptionsInt("-mat_strumpack_hss_min_sep_size","Minimum size of separator for HSS compression","None",hssminsize,&hssminsize,&set);CHKERRQ(ierr);
579291d0ed5SPieter Ghysels   if (set) PetscStackCall("STRUMPACK_set_HSS_min_sep_size",STRUMPACK_set_HSS_min_sep_size(*S,(int)hssminsize));
580575704cbSPieter Ghysels 
581a36bf211SPieter Ghysels   PetscStackCall("STRUMPACK_HSS_max_rank",max_rank = (PetscInt)STRUMPACK_HSS_max_rank(*S));
58234f43fa5SPieter Ghysels   ierr = PetscOptionsInt("-mat_strumpack_max_rank","Maximum rank in HSS approximation","None",max_rank,&max_rank,&set);CHKERRQ(ierr);
583a36bf211SPieter Ghysels   if (set) PetscStackCall("STRUMPACK_set_HSS_max_rank",STRUMPACK_set_HSS_max_rank(*S,(int)max_rank));
584a36bf211SPieter Ghysels 
585a36bf211SPieter Ghysels   PetscStackCall("STRUMPACK_HSS_leaf_size",leaf_size = (PetscInt)STRUMPACK_HSS_leaf_size(*S));
586a36bf211SPieter Ghysels   ierr = PetscOptionsInt("-mat_strumpack_leaf_size","Size of diagonal blocks in HSS approximation","None",leaf_size,&leaf_size,&set);CHKERRQ(ierr);
587a36bf211SPieter Ghysels   if (set) PetscStackCall("STRUMPACK_set_HSS_leaf_size",STRUMPACK_set_HSS_leaf_size(*S,(int)leaf_size));
58834f43fa5SPieter Ghysels 
58934f43fa5SPieter Ghysels   const char *const STRUMPACKNDTypes[] = {"NATURAL","METIS","PARMETIS","SCOTCH","PTSCOTCH","RCM","STRUMPACKNDTypes","",0};
59034f43fa5SPieter Ghysels   PetscStackCall("STRUMPACK_reordering_method",ndcurrent = STRUMPACK_reordering_method(*S));
59134f43fa5SPieter Ghysels   PetscOptionsEnum("-mat_strumpack_reordering","Sparsity reducing matrix reordering","None",STRUMPACKNDTypes,(PetscEnum)ndcurrent,(PetscEnum*)&ndvalue,&set);CHKERRQ(ierr);
59234f43fa5SPieter Ghysels   if (set) PetscStackCall("STRUMPACK_set_reordering_method",STRUMPACK_set_reordering_method(*S,ndvalue));
593575704cbSPieter Ghysels 
59488382b05SPieter Ghysels   if (ftype == MAT_FACTOR_ILU) {
59588382b05SPieter Ghysels     /* When enabling HSS compression, the STRUMPACK solver becomes an incomplete              */
59688382b05SPieter Ghysels     /* (or approximate) LU factorization.                                                     */
597291d0ed5SPieter Ghysels     PetscStackCall("STRUMPACK_enable_HSS",STRUMPACK_enable_HSS(*S));
59888382b05SPieter Ghysels   }
599575704cbSPieter Ghysels 
60034f43fa5SPieter Ghysels   /* Disable the outer iterative solver from STRUMPACK.                                       */
60134f43fa5SPieter Ghysels   /* When STRUMPACK is used as a direct solver, it will by default do iterative refinement.   */
60234f43fa5SPieter Ghysels   /* When STRUMPACK is used as an approximate factorization preconditioner (by enabling       */
60334f43fa5SPieter Ghysels   /* low-rank compression), it will use it's own preconditioned GMRES. Here we can disable    */
60434f43fa5SPieter Ghysels   /* the outer iterative solver, as PETSc uses STRUMPACK from within a KSP.                   */
60534f43fa5SPieter Ghysels   PetscStackCall("STRUMPACK_set_Krylov_solver", STRUMPACK_set_Krylov_solver(*S, STRUMPACK_DIRECT));
60634f43fa5SPieter Ghysels 
60734f43fa5SPieter Ghysels   const char *const SolverTypes[] = {"AUTO","NONE","REFINE","PREC_GMRES","GMRES","PREC_BICGSTAB","BICGSTAB","SolverTypes","",0};
60834f43fa5SPieter Ghysels   PetscStackCall("STRUMPACK_Krylov_solver",itcurrent = STRUMPACK_Krylov_solver(*S));
60934f43fa5SPieter Ghysels   PetscOptionsEnum("-mat_strumpack_iterative_solver","Select iterative solver from STRUMPACK","None",SolverTypes,(PetscEnum)itcurrent,(PetscEnum*)&itsolver,&set);CHKERRQ(ierr);
61034f43fa5SPieter Ghysels   if (set) PetscStackCall("STRUMPACK_set_Krylov_solver",STRUMPACK_set_Krylov_solver(*S,itsolver));
61134f43fa5SPieter Ghysels 
61234f43fa5SPieter Ghysels   PetscOptionsEnd();
61355c022e5SPieter Ghysels 
6147d6ea485SPieter Ghysels   *F = B;
6157d6ea485SPieter Ghysels   PetscFunctionReturn(0);
6167d6ea485SPieter Ghysels }
6177d6ea485SPieter Ghysels 
618291d0ed5SPieter Ghysels #undef __FUNCT__
619291d0ed5SPieter Ghysels #define __FUNCT__ "MatSolverPackageRegister_STRUMPACK"
6207d6ea485SPieter Ghysels PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_STRUMPACK(void)
6217d6ea485SPieter Ghysels {
6227d6ea485SPieter Ghysels   PetscErrorCode ierr;
6237d6ea485SPieter Ghysels 
6247d6ea485SPieter Ghysels   PetscFunctionBegin;
6257d6ea485SPieter Ghysels   ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr);
6267d6ea485SPieter Ghysels   ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr);
627575704cbSPieter Ghysels   ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_ILU,MatGetFactor_aij_strumpack);CHKERRQ(ierr);
628575704cbSPieter Ghysels   ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_ILU,MatGetFactor_aij_strumpack);CHKERRQ(ierr);
6297d6ea485SPieter Ghysels   PetscFunctionReturn(0);
6307d6ea485SPieter Ghysels }
631