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