xref: /petsc/src/mat/impls/aij/mpi/strumpack/strumpack.c (revision f155c23239e6d3f7c7ec79ff00b4f28519d0ce99)
1 #include <../src/mat/impls/aij/seq/aij.h>            /*I "petscmat.h" I*/
2 #include <../src/mat/impls/aij/mpi/mpiaij.h>
3 #include <StrumpackSparseSolver.h>
4 
5 static PetscErrorCode MatGetDiagonal_STRUMPACK(Mat A,Vec v)
6 {
7   PetscFunctionBegin;
8   SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type: STRUMPACK factor");
9   PetscFunctionReturn(0);
10 }
11 
12 static PetscErrorCode MatDestroy_STRUMPACK(Mat A)
13 {
14   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)A->spptr;
15   PetscErrorCode         ierr;
16   PetscBool              flg;
17 
18   PetscFunctionBegin;
19   /* Deallocate STRUMPACK storage */
20   PetscStackCall("STRUMPACK_destroy",STRUMPACK_destroy(S));
21   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
22   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr);
23   if (flg) {
24     ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
25   } else {
26     ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr);
27   }
28 
29   /* clear composed functions */
30   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverType_C",NULL);CHKERRQ(ierr);
31   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetReordering_C",NULL);CHKERRQ(ierr);
32   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetColPerm_C",NULL);CHKERRQ(ierr);
33   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSRelTol_C",NULL);CHKERRQ(ierr);
34   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSAbsTol_C",NULL);CHKERRQ(ierr);
35   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSMaxRank_C",NULL);CHKERRQ(ierr);
36   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSLeafSize_C",NULL);CHKERRQ(ierr);
37   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSMinSepSize_C",NULL);CHKERRQ(ierr);
38 
39   PetscFunctionReturn(0);
40 }
41 
42 static PetscErrorCode MatSTRUMPACKSetReordering_STRUMPACK(Mat F,MatSTRUMPACKReordering reordering)
43 {
44   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
45 
46   PetscFunctionBegin;
47   PetscStackCall("STRUMPACK_reordering_method",STRUMPACK_set_reordering_method(*S,(STRUMPACK_REORDERING_STRATEGY)reordering));
48   PetscFunctionReturn(0);
49 }
50 
51 /*@
52   MatSTRUMPACKSetReordering - Set STRUMPACK fill-reducing reordering
53 
54    Input Parameters:
55 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
56 -  reordering - the code to be used to find the fill-reducing reordering
57       Possible values: NATURAL=0 METIS=1 PARMETIS=2 SCOTCH=3 PTSCOTCH=4 RCM=5
58 
59   Options Database:
60 .   -mat_strumpack_reordering <METIS>  - Sparsity reducing matrix reordering (choose one of) NATURAL METIS PARMETIS SCOTCH PTSCOTCH RCM (None)
61 
62    Level: beginner
63 
64    References:
65 .  * - STRUMPACK manual
66 
67 .seealso: MatGetFactor()
68 @*/
69 PetscErrorCode MatSTRUMPACKSetReordering(Mat F,MatSTRUMPACKReordering reordering)
70 {
71   PetscErrorCode ierr;
72 
73   PetscFunctionBegin;
74   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
75   PetscValidLogicalCollectiveEnum(F,reordering,2);
76   ierr = PetscTryMethod(F,"MatSTRUMPACKSetReordering_C",(Mat,MatSTRUMPACKReordering),(F,reordering));CHKERRQ(ierr);
77   PetscFunctionReturn(0);
78 }
79 
80 static PetscErrorCode MatSTRUMPACKSetColPerm_STRUMPACK(Mat F,PetscBool cperm)
81 {
82   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
83 
84   PetscFunctionBegin;
85   PetscStackCall("STRUMPACK_set_mc64job",STRUMPACK_set_mc64job(*S,cperm ? 5 : 0));
86   PetscFunctionReturn(0);
87 }
88 
89 /*@
90   MatSTRUMPACKSetColPerm - Set whether STRUMPACK should try to permute the columns of the matrix in order to get a nonzero diagonal
91 
92    Logically Collective on Mat
93 
94    Input Parameters:
95 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
96 -  cperm - PETSC_TRUE to permute (internally) the columns of the matrix
97 
98   Options Database:
99 .   -mat_strumpack_colperm <cperm> - true to use the permutation
100 
101    Level: beginner
102 
103    References:
104 .  * - STRUMPACK manual
105 
106 .seealso: MatGetFactor()
107 @*/
108 PetscErrorCode MatSTRUMPACKSetColPerm(Mat F,PetscBool cperm)
109 {
110   PetscErrorCode ierr;
111 
112   PetscFunctionBegin;
113   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
114   PetscValidLogicalCollectiveBool(F,cperm,2);
115   ierr = PetscTryMethod(F,"MatSTRUMPACKSetColPerm_C",(Mat,PetscBool),(F,cperm));CHKERRQ(ierr);
116   PetscFunctionReturn(0);
117 }
118 
119 static PetscErrorCode MatSTRUMPACKSetHSSRelTol_STRUMPACK(Mat F,PetscReal rtol)
120 {
121   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
122 
123   PetscFunctionBegin;
124   PetscStackCall("STRUMPACK_set_HSS_rel_tol", STRUMPACK_set_HSS_rel_tol(*S,rtol));
125   PetscFunctionReturn(0);
126 }
127 
128 /*@
129   MatSTRUMPACKSetHSSRelTol - Set STRUMPACK relative tolerance for HSS compression
130 
131   Logically Collective on Mat
132 
133    Input Parameters:
134 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
135 -  rtol - relative compression tolerance
136 
137   Options Database:
138 .   -mat_strumpack_hss_rel_tol <1e-2>         - Relative compression tolerance (None)
139 
140    Level: beginner
141 
142    References:
143 .  * - STRUMPACK manual
144 
145 .seealso: MatGetFactor()
146 @*/
147 PetscErrorCode MatSTRUMPACKSetHSSRelTol(Mat F,PetscReal rtol)
148 {
149   PetscErrorCode ierr;
150 
151   PetscFunctionBegin;
152   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
153   PetscValidLogicalCollectiveReal(F,rtol,2);
154   ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSRelTol_C",(Mat,PetscReal),(F,rtol));CHKERRQ(ierr);
155   PetscFunctionReturn(0);
156 }
157 
158 static PetscErrorCode MatSTRUMPACKSetHSSAbsTol_STRUMPACK(Mat F,PetscReal atol)
159 {
160   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
161 
162   PetscFunctionBegin;
163   PetscStackCall("STRUMPACK_set_HSS_abs_tol", STRUMPACK_set_HSS_abs_tol(*S,atol));
164   PetscFunctionReturn(0);
165 }
166 
167 /*@
168   MatSTRUMPACKSetHSSAbsTol - Set STRUMPACK absolute tolerance for HSS compression
169 
170    Logically Collective on Mat
171 
172    Input Parameters:
173 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
174 -  atol - absolute compression tolerance
175 
176   Options Database:
177 .   -mat_strumpack_hss_abs_tol <1e-8>         - Absolute compression tolerance (None)
178 
179    Level: beginner
180 
181    References:
182 .  * - STRUMPACK manual
183 
184 .seealso: MatGetFactor()
185 @*/
186 PetscErrorCode MatSTRUMPACKSetHSSAbsTol(Mat F,PetscReal atol)
187 {
188   PetscErrorCode ierr;
189 
190   PetscFunctionBegin;
191   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
192   PetscValidLogicalCollectiveReal(F,atol,2);
193   ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSAbsTol_C",(Mat,PetscReal),(F,atol));CHKERRQ(ierr);
194   PetscFunctionReturn(0);
195 }
196 
197 static PetscErrorCode MatSTRUMPACKSetHSSMaxRank_STRUMPACK(Mat F,PetscInt hssmaxrank)
198 {
199   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
200 
201   PetscFunctionBegin;
202   PetscStackCall("STRUMPACK_set_HSS_max_rank", STRUMPACK_set_HSS_max_rank(*S,hssmaxrank));
203   PetscFunctionReturn(0);
204 }
205 
206 /*@
207   MatSTRUMPACKSetHSSMaxRank - Set STRUMPACK maximum HSS rank
208 
209    Logically Collective on Mat
210 
211    Input Parameters:
212 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
213 -  hssmaxrank - maximum rank used in low-rank approximation
214 
215   Options Database:
216 .   -mat_strumpack_max_rank    - Maximum rank in HSS compression, when using pctype ilu (None)
217 
218    Level: beginner
219 
220    References:
221 .  * - STRUMPACK manual
222 
223 .seealso: MatGetFactor()
224 @*/
225 PetscErrorCode MatSTRUMPACKSetHSSMaxRank(Mat F,PetscInt hssmaxrank)
226 {
227   PetscErrorCode ierr;
228 
229   PetscFunctionBegin;
230   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
231   PetscValidLogicalCollectiveInt(F,hssmaxrank,2);
232   ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSMaxRank_C",(Mat,PetscInt),(F,hssmaxrank));CHKERRQ(ierr);
233   PetscFunctionReturn(0);
234 }
235 
236 static PetscErrorCode MatSTRUMPACKSetHSSLeafSize_STRUMPACK(Mat F,PetscInt leaf_size)
237 {
238   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
239 
240   PetscFunctionBegin;
241   PetscStackCall("STRUMPACK_set_HSS_leaf_size", STRUMPACK_set_HSS_leaf_size(*S,leaf_size));
242   PetscFunctionReturn(0);
243 }
244 
245 /*@
246   MatSTRUMPACKSetHSSLeafSize - Set STRUMPACK HSS leaf size
247 
248    Logically Collective on Mat
249 
250    Input Parameters:
251 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
252 -  leaf_size - Size of diagonal blocks in HSS approximation
253 
254   Options Database:
255 .   -mat_strumpack_leaf_size    - Size of diagonal blocks in HSS approximation, when using pctype ilu (None)
256 
257    Level: beginner
258 
259    References:
260 .  * - STRUMPACK manual
261 
262 .seealso: MatGetFactor()
263 @*/
264 PetscErrorCode MatSTRUMPACKSetHSSLeafSize(Mat F,PetscInt leaf_size)
265 {
266   PetscErrorCode ierr;
267 
268   PetscFunctionBegin;
269   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
270   PetscValidLogicalCollectiveInt(F,leaf_size,2);
271   ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSLeafSize_C",(Mat,PetscInt),(F,leaf_size));CHKERRQ(ierr);
272   PetscFunctionReturn(0);
273 }
274 
275 static PetscErrorCode MatSTRUMPACKSetHSSMinSepSize_STRUMPACK(Mat F,PetscInt hssminsize)
276 {
277   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
278 
279   PetscFunctionBegin;
280   PetscStackCall("STRUMPACK_set_HSS_min_sep_size", STRUMPACK_set_HSS_min_sep_size(*S,hssminsize));
281   PetscFunctionReturn(0);
282 }
283 
284 /*@
285   MatSTRUMPACKSetHSSMinSepSize - Set STRUMPACK minimum separator size for low-rank approximation
286 
287    Logically Collective on Mat
288 
289    Input Parameters:
290 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
291 -  hssminsize - minimum dense matrix size for low-rank approximation
292 
293   Options Database:
294 .   -mat_strumpack_hss_min_sep_size <hssminsize> - set the minimum separator size
295 
296    Level: beginner
297 
298    References:
299 .  * - STRUMPACK manual
300 
301 .seealso: MatGetFactor()
302 @*/
303 PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat F,PetscInt hssminsize)
304 {
305   PetscErrorCode ierr;
306 
307   PetscFunctionBegin;
308   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
309   PetscValidLogicalCollectiveInt(F,hssminsize,2);
310   ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSMinSepSize_C",(Mat,PetscInt),(F,hssminsize));CHKERRQ(ierr);
311   PetscFunctionReturn(0);
312 }
313 
314 static PetscErrorCode MatSolve_STRUMPACK(Mat A,Vec b_mpi,Vec x)
315 {
316   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)A->spptr;
317   STRUMPACK_RETURN_CODE  sp_err;
318   PetscErrorCode         ierr;
319   const PetscScalar      *bptr;
320   PetscScalar            *xptr;
321 
322   PetscFunctionBegin;
323   ierr = VecGetArray(x,&xptr);CHKERRQ(ierr);
324   ierr = VecGetArrayRead(b_mpi,&bptr);CHKERRQ(ierr);
325 
326   PetscStackCall("STRUMPACK_solve",sp_err = STRUMPACK_solve(*S,(PetscScalar*)bptr,xptr,0));
327   switch (sp_err) {
328   case STRUMPACK_SUCCESS: break;
329   case STRUMPACK_MATRIX_NOT_SET:   { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); break; }
330   case STRUMPACK_REORDERING_ERROR: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); break; }
331   default:                           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: solve failed");
332   }
333   ierr = VecRestoreArray(x,&xptr);CHKERRQ(ierr);
334   ierr = VecRestoreArrayRead(b_mpi,&bptr);CHKERRQ(ierr);
335   PetscFunctionReturn(0);
336 }
337 
338 static PetscErrorCode MatMatSolve_STRUMPACK(Mat A,Mat B_mpi,Mat X)
339 {
340   PetscErrorCode   ierr;
341   PetscBool        flg;
342 
343   PetscFunctionBegin;
344   ierr = PetscObjectTypeCompareAny((PetscObject)B_mpi,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
345   PetscCheckFalse(!flg,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
346   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
347   PetscCheckFalse(!flg,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
348   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_STRUMPACK() is not implemented yet");
349   PetscFunctionReturn(0);
350 }
351 
352 static PetscErrorCode MatView_Info_STRUMPACK(Mat A,PetscViewer viewer)
353 {
354   PetscErrorCode  ierr;
355 
356   PetscFunctionBegin;
357   /* check if matrix is strumpack type */
358   if (A->ops->solve != MatSolve_STRUMPACK) PetscFunctionReturn(0);
359   ierr = PetscViewerASCIIPrintf(viewer,"STRUMPACK sparse solver!\n");CHKERRQ(ierr);
360   PetscFunctionReturn(0);
361 }
362 
363 static PetscErrorCode MatView_STRUMPACK(Mat A,PetscViewer viewer)
364 {
365   PetscErrorCode    ierr;
366   PetscBool         iascii;
367   PetscViewerFormat format;
368 
369   PetscFunctionBegin;
370   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
371   if (iascii) {
372     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
373     if (format == PETSC_VIEWER_ASCII_INFO) {
374       ierr = MatView_Info_STRUMPACK(A,viewer);CHKERRQ(ierr);
375     }
376   }
377   PetscFunctionReturn(0);
378 }
379 
380 static PetscErrorCode MatLUFactorNumeric_STRUMPACK(Mat F,Mat A,const MatFactorInfo *info)
381 {
382   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
383   STRUMPACK_RETURN_CODE  sp_err;
384   Mat_SeqAIJ             *A_d,*A_o;
385   Mat_MPIAIJ             *mat;
386   PetscErrorCode         ierr;
387   PetscInt               M=A->rmap->N,m=A->rmap->n;
388   PetscBool              flg;
389 
390   PetscFunctionBegin;
391   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
392   if (flg) { /* A is MATMPIAIJ */
393     mat = (Mat_MPIAIJ*)A->data;
394     A_d = (Mat_SeqAIJ*)(mat->A)->data;
395     A_o = (Mat_SeqAIJ*)(mat->B)->data;
396     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));
397   } else { /* A is MATSEQAIJ */
398     A_d = (Mat_SeqAIJ*)A->data;
399     PetscStackCall("STRUMPACK_set_csr_matrix",STRUMPACK_set_csr_matrix(*S,&M,A_d->i,A_d->j,A_d->a,0));
400   }
401 
402   /* Reorder and Factor the matrix. */
403   /* TODO figure out how to avoid reorder if the matrix values changed, but the pattern remains the same. */
404   PetscStackCall("STRUMPACK_reorder",sp_err = STRUMPACK_reorder(*S));
405   PetscStackCall("STRUMPACK_factor",sp_err = STRUMPACK_factor(*S));
406   switch (sp_err) {
407   case STRUMPACK_SUCCESS: break;
408   case STRUMPACK_MATRIX_NOT_SET:   { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); break; }
409   case STRUMPACK_REORDERING_ERROR: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); break; }
410   default:                           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: factorization failed");
411   }
412   F->assembled    = PETSC_TRUE;
413   F->preallocated = PETSC_TRUE;
414   PetscFunctionReturn(0);
415 }
416 
417 static PetscErrorCode MatLUFactorSymbolic_STRUMPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
418 {
419   PetscFunctionBegin;
420   F->ops->lufactornumeric = MatLUFactorNumeric_STRUMPACK;
421   F->ops->solve           = MatSolve_STRUMPACK;
422   F->ops->matsolve        = MatMatSolve_STRUMPACK;
423   PetscFunctionReturn(0);
424 }
425 
426 static PetscErrorCode MatFactorGetSolverType_aij_strumpack(Mat A,MatSolverType *type)
427 {
428   PetscFunctionBegin;
429   *type = MATSOLVERSTRUMPACK;
430   PetscFunctionReturn(0);
431 }
432 
433 /*MC
434   MATSOLVERSSTRUMPACK = "strumpack" - A solver package providing a direct sparse solver (PCLU)
435   and a preconditioner (PCILU) using low-rank compression via the external package STRUMPACK.
436 
437   Consult the STRUMPACK-sparse manual for more info.
438 
439   Use
440      ./configure --download-strumpack
441   to have PETSc installed with STRUMPACK
442 
443   Use
444     -pc_type lu -pc_factor_mat_solver_type strumpack
445   to use this as an exact (direct) solver, use
446     -pc_type ilu -pc_factor_mat_solver_type strumpack
447   to enable low-rank compression (i.e, use as a preconditioner).
448 
449   Works with AIJ matrices
450 
451   Options Database Keys:
452 + -mat_strumpack_verbose                    - verbose info
453 . -mat_strumpack_hss_rel_tol <1e-2>         - Relative compression tolerance (None)
454 . -mat_strumpack_hss_abs_tol <1e-8>         - Absolute compression tolerance (None)
455 . -mat_strumpack_colperm <TRUE>             - Permute matrix to make diagonal nonzeros (None)
456 . -mat_strumpack_hss_min_sep_size <256>     - Minimum size of separator for HSS compression (None)
457 . -mat_strumpack_max_rank                   - Maximum rank in HSS compression, when using pctype ilu (None)
458 . -mat_strumpack_leaf_size                  - Size of diagonal blocks in HSS approximation, when using pctype ilu (None)
459 . -mat_strumpack_reordering <METIS>         - Sparsity reducing matrix reordering (choose one of) NATURAL METIS PARMETIS SCOTCH PTSCOTCH RCM (None)
460 - -mat_strumpack_iterative_solver <DIRECT>  - Select iterative solver from STRUMPACK (choose one of) AUTO DIRECT REFINE PREC_GMRES GMRES PREC_BICGSTAB BICGSTAB (None)
461 
462  Level: beginner
463 
464 .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, PCFactorSetMatSolverType(), MatSolverType
465 M*/
466 static PetscErrorCode MatGetFactor_aij_strumpack(Mat A,MatFactorType ftype,Mat *F)
467 {
468   Mat                           B;
469   PetscErrorCode                ierr;
470   PetscInt                      M=A->rmap->N,N=A->cmap->N;
471   PetscBool                     verb,flg,set;
472   PetscReal                     ctol;
473   PetscInt                      hssminsize,max_rank,leaf_size;
474   STRUMPACK_SparseSolver        *S;
475   STRUMPACK_INTERFACE           iface;
476   STRUMPACK_REORDERING_STRATEGY ndcurrent,ndvalue;
477   STRUMPACK_KRYLOV_SOLVER       itcurrent,itsolver;
478   const STRUMPACK_PRECISION     table[2][2][2] =
479     {{{STRUMPACK_FLOATCOMPLEX_64, STRUMPACK_DOUBLECOMPLEX_64},
480       {STRUMPACK_FLOAT_64,        STRUMPACK_DOUBLE_64}},
481      {{STRUMPACK_FLOATCOMPLEX,    STRUMPACK_DOUBLECOMPLEX},
482       {STRUMPACK_FLOAT,           STRUMPACK_DOUBLE}}};
483   const STRUMPACK_PRECISION     prec = table[(sizeof(PetscInt)==8)?0:1][(PETSC_SCALAR==PETSC_COMPLEX)?0:1][(PETSC_REAL==PETSC_FLOAT)?0:1];
484   const char *const             STRUMPACKNDTypes[] = {"NATURAL","METIS","PARMETIS","SCOTCH","PTSCOTCH","RCM","STRUMPACKNDTypes","",0};
485   const char *const             SolverTypes[] = {"AUTO","NONE","REFINE","PREC_GMRES","GMRES","PREC_BICGSTAB","BICGSTAB","SolverTypes","",0};
486 
487   PetscFunctionBegin;
488   /* Create the factorization matrix */
489   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
490   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,M,N);CHKERRQ(ierr);
491   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
492   ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr);
493   ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr);
494   B->trivialsymbolic = PETSC_TRUE;
495   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) {
496     B->ops->lufactorsymbolic  = MatLUFactorSymbolic_STRUMPACK;
497     B->ops->ilufactorsymbolic = MatLUFactorSymbolic_STRUMPACK;
498   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
499   B->ops->view        = MatView_STRUMPACK;
500   B->ops->destroy     = MatDestroy_STRUMPACK;
501   B->ops->getdiagonal = MatGetDiagonal_STRUMPACK;
502   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_aij_strumpack);CHKERRQ(ierr);
503   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetReordering_C",MatSTRUMPACKSetReordering_STRUMPACK);CHKERRQ(ierr);
504   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetColPerm_C",MatSTRUMPACKSetColPerm_STRUMPACK);CHKERRQ(ierr);
505   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSRelTol_C",MatSTRUMPACKSetHSSRelTol_STRUMPACK);CHKERRQ(ierr);
506   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSAbsTol_C",MatSTRUMPACKSetHSSAbsTol_STRUMPACK);CHKERRQ(ierr);
507   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSMaxRank_C",MatSTRUMPACKSetHSSMaxRank_STRUMPACK);CHKERRQ(ierr);
508   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSLeafSize_C",MatSTRUMPACKSetHSSLeafSize_STRUMPACK);CHKERRQ(ierr);
509   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSMinSepSize_C",MatSTRUMPACKSetHSSMinSepSize_STRUMPACK);CHKERRQ(ierr);
510   B->factortype = ftype;
511   ierr     = PetscNewLog(B,&S);CHKERRQ(ierr);
512   B->spptr = S;
513 
514   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr);
515   iface = flg ? STRUMPACK_MT : STRUMPACK_MPI_DIST;
516 
517   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"STRUMPACK Options","Mat");CHKERRQ(ierr);
518 
519   verb = PetscLogPrintInfo ? PETSC_TRUE : PETSC_FALSE;
520   ierr = PetscOptionsBool("-mat_strumpack_verbose","Print STRUMPACK information","None",verb,&verb,NULL);CHKERRQ(ierr);
521 
522   PetscStackCall("STRUMPACK_init",STRUMPACK_init(S,PetscObjectComm((PetscObject)A),prec,iface,0,NULL,verb));
523 
524   PetscStackCall("STRUMPACK_HSS_rel_tol",ctol = (PetscReal)STRUMPACK_HSS_rel_tol(*S));
525   ierr = PetscOptionsReal("-mat_strumpack_hss_rel_tol","Relative compression tolerance","None",ctol,&ctol,&set);CHKERRQ(ierr);
526   if (set) PetscStackCall("STRUMPACK_set_HSS_rel_tol",STRUMPACK_set_HSS_rel_tol(*S,(double)ctol));
527 
528   PetscStackCall("STRUMPACK_HSS_abs_tol",ctol = (PetscReal)STRUMPACK_HSS_abs_tol(*S));
529   ierr = PetscOptionsReal("-mat_strumpack_hss_abs_tol","Absolute compression tolerance","None",ctol,&ctol,&set);CHKERRQ(ierr);
530   if (set) PetscStackCall("STRUMPACK_set_HSS_abs_tol",STRUMPACK_set_HSS_abs_tol(*S,(double)ctol));
531 
532   PetscStackCall("STRUMPACK_mc64job",flg = (STRUMPACK_mc64job(*S) == 0) ? PETSC_FALSE : PETSC_TRUE);
533   ierr = PetscOptionsBool("-mat_strumpack_colperm","Find a col perm to get nonzero diagonal","None",flg,&flg,&set);CHKERRQ(ierr);
534   if (set) PetscStackCall("STRUMPACK_set_mc64job",STRUMPACK_set_mc64job(*S,flg ? 5 : 0));
535 
536   PetscStackCall("STRUMPACK_HSS_min_sep_size",hssminsize = (PetscInt)STRUMPACK_HSS_min_sep_size(*S));
537   ierr = PetscOptionsInt("-mat_strumpack_hss_min_sep_size","Minimum size of separator for HSS compression","None",hssminsize,&hssminsize,&set);CHKERRQ(ierr);
538   if (set) PetscStackCall("STRUMPACK_set_HSS_min_sep_size",STRUMPACK_set_HSS_min_sep_size(*S,(int)hssminsize));
539 
540   PetscStackCall("STRUMPACK_HSS_max_rank",max_rank = (PetscInt)STRUMPACK_HSS_max_rank(*S));
541   ierr = PetscOptionsInt("-mat_strumpack_max_rank","Maximum rank in HSS approximation","None",max_rank,&max_rank,&set);CHKERRQ(ierr);
542   if (set) PetscStackCall("STRUMPACK_set_HSS_max_rank",STRUMPACK_set_HSS_max_rank(*S,(int)max_rank));
543 
544   PetscStackCall("STRUMPACK_HSS_leaf_size",leaf_size = (PetscInt)STRUMPACK_HSS_leaf_size(*S));
545   ierr = PetscOptionsInt("-mat_strumpack_leaf_size","Size of diagonal blocks in HSS approximation","None",leaf_size,&leaf_size,&set);CHKERRQ(ierr);
546   if (set) PetscStackCall("STRUMPACK_set_HSS_leaf_size",STRUMPACK_set_HSS_leaf_size(*S,(int)leaf_size));
547 
548   PetscStackCall("STRUMPACK_reordering_method",ndcurrent = STRUMPACK_reordering_method(*S));
549   ierr = PetscOptionsEnum("-mat_strumpack_reordering","Sparsity reducing matrix reordering","None",STRUMPACKNDTypes,(PetscEnum)ndcurrent,(PetscEnum*)&ndvalue,&set);CHKERRQ(ierr);
550   if (set) PetscStackCall("STRUMPACK_set_reordering_method",STRUMPACK_set_reordering_method(*S,ndvalue));
551 
552   if (ftype == MAT_FACTOR_ILU) {
553     /* When enabling HSS compression, the STRUMPACK solver becomes an incomplete              */
554     /* (or approximate) LU factorization.                                                     */
555     PetscStackCall("STRUMPACK_enable_HSS",STRUMPACK_enable_HSS(*S));
556   }
557 
558   /* Disable the outer iterative solver from STRUMPACK.                                       */
559   /* When STRUMPACK is used as a direct solver, it will by default do iterative refinement.   */
560   /* When STRUMPACK is used as an approximate factorization preconditioner (by enabling       */
561   /* low-rank compression), it will use it's own preconditioned GMRES. Here we can disable    */
562   /* the outer iterative solver, as PETSc uses STRUMPACK from within a KSP.                   */
563   PetscStackCall("STRUMPACK_set_Krylov_solver", STRUMPACK_set_Krylov_solver(*S, STRUMPACK_DIRECT));
564 
565   PetscStackCall("STRUMPACK_Krylov_solver",itcurrent = STRUMPACK_Krylov_solver(*S));
566   ierr = PetscOptionsEnum("-mat_strumpack_iterative_solver","Select iterative solver from STRUMPACK","None",SolverTypes,(PetscEnum)itcurrent,(PetscEnum*)&itsolver,&set);CHKERRQ(ierr);
567   if (set) PetscStackCall("STRUMPACK_set_Krylov_solver",STRUMPACK_set_Krylov_solver(*S,itsolver));
568 
569   PetscOptionsEnd();
570 
571   *F = B;
572   PetscFunctionReturn(0);
573 }
574 
575 PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_STRUMPACK(void)
576 {
577   PetscErrorCode ierr;
578 
579   PetscFunctionBegin;
580   ierr = MatSolverTypeRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr);
581   ierr = MatSolverTypeRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr);
582   ierr = MatSolverTypeRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_ILU,MatGetFactor_aij_strumpack);CHKERRQ(ierr);
583   ierr = MatSolverTypeRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_ILU,MatGetFactor_aij_strumpack);CHKERRQ(ierr);
584   PetscFunctionReturn(0);
585 }
586