xref: /petsc/src/mat/impls/aij/mpi/strumpack/strumpack.c (revision 34f43fa5948b3ce7830866b2340111e5aa7acbcc)
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,"MatSTRUMPACKSetHSSMinFrontSize_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_max_rank", STRUMPACK_set_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__ "MatSTRUMPACKSetHSSMinFrontSize_STRUMPACK"
258 static PetscErrorCode MatSTRUMPACKSetHSSMinFrontSize_STRUMPACK(Mat F,PetscInt hssminsize)
259 {
260   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
261 
262   PetscFunctionBegin;
263   PetscStackCall("STRUMPACK_set_HSS_min_front_size", STRUMPACK_set_HSS_min_front_size(*S,hssminsize));
264   PetscFunctionReturn(0);
265 }
266 #undef __FUNCT__
267 #define __FUNCT__ "MatSTRUMPACKSetHSSMinFrontSize"
268 /*@
269   MatSTRUMPACKSetHSSMinFrontSize - Set STRUMPACK minimum dense matrix size for low-rank approximation
270    Logically Collective on Mat
271 
272    Input Parameters:
273 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
274 -  hssminsize - minimum dense matrix size for low-rank approximation
275 
276   Options Database:
277 .   -mat_strumpack_hss_min_front_size <hssminsize>
278 
279    Level: beginner
280 
281    References:
282 .      STRUMPACK manual
283 
284 .seealso: MatGetFactor()
285 @*/
286 PetscErrorCode MatSTRUMPACKSetHSSMinFrontSize(Mat F,PetscInt hssminsize)
287 {
288   PetscErrorCode ierr;
289 
290   PetscFunctionBegin;
291   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
292   PetscValidLogicalCollectiveInt(F,hssminsize,2);
293   ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSMinFrontSize_C",(Mat,PetscInt),(F,hssminsize));CHKERRQ(ierr);
294   PetscFunctionReturn(0);
295 }
296 
297 #undef __FUNCT__
298 #define __FUNCT__ "MatSTRUMPACKSetHSSMinSepSize_STRUMPACK"
299 static PetscErrorCode MatSTRUMPACKSetHSSMinSepSize_STRUMPACK(Mat F,PetscInt hssminsize)
300 {
301   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
302 
303   PetscFunctionBegin;
304   PetscStackCall("STRUMPACK_set_HSS_min_sep_size", STRUMPACK_set_HSS_min_sep_size(*S,hssminsize));
305   PetscFunctionReturn(0);
306 }
307 #undef __FUNCT__
308 #define __FUNCT__ "MatSTRUMPACKSetHSSMinSepSize"
309 /*@
310   MatSTRUMPACKSetHSSMinSepSize - Set STRUMPACK minimum separator size for low-rank approximation
311    Logically Collective on Mat
312 
313    Input Parameters:
314 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
315 -  hssminsize - minimum dense matrix size for low-rank approximation
316 
317   Options Database:
318 .   -mat_strumpack_hss_min_sep_size <hssminsize>
319 
320    Level: beginner
321 
322    References:
323 .      STRUMPACK manual
324 
325 .seealso: MatGetFactor()
326 @*/
327 PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat F,PetscInt hssminsize)
328 {
329   PetscErrorCode ierr;
330 
331   PetscFunctionBegin;
332   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
333   PetscValidLogicalCollectiveInt(F,hssminsize,2);
334   ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSMinSepSize_C",(Mat,PetscInt),(F,hssminsize));CHKERRQ(ierr);
335   PetscFunctionReturn(0);
336 }
337 
338 
339 #undef __FUNCT__
340 #define __FUNCT__ "MatSolve_STRUMPACK"
341 static PetscErrorCode MatSolve_STRUMPACK(Mat A,Vec b_mpi,Vec x)
342 {
343   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)A->spptr;
344   STRUMPACK_RETURN_CODE  sp_err;
345   PetscErrorCode         ierr;
346   const PetscScalar      *bptr;
347   PetscScalar            *xptr;
348 
349   PetscFunctionBegin;
350   ierr = VecGetArray(x,&xptr);CHKERRQ(ierr);
351   ierr = VecGetArrayRead(b_mpi,&bptr);CHKERRQ(ierr);
352 
353   PetscStackCall("STRUMPACK_solve",sp_err = STRUMPACK_solve(*S,(PetscScalar*)bptr,xptr,0));
354   switch (sp_err) {
355   case STRUMPACK_SUCCESS: break;
356   case STRUMPACK_MATRIX_NOT_SET:   { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); break; }
357   case STRUMPACK_REORDERING_ERROR: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); break; }
358   default:                           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: solve failed");
359   }
360   ierr = VecRestoreArray(x,&xptr);CHKERRQ(ierr);
361   ierr = VecRestoreArrayRead(b_mpi,&bptr);CHKERRQ(ierr);
362   PetscFunctionReturn(0);
363 }
364 
365 #undef __FUNCT__
366 #define __FUNCT__ "MatMatSolve_STRUMPACK"
367 static PetscErrorCode MatMatSolve_STRUMPACK(Mat A,Mat B_mpi,Mat X)
368 {
369   PetscErrorCode   ierr;
370   PetscBool        flg;
371 
372   PetscFunctionBegin;
373   ierr = PetscObjectTypeCompareAny((PetscObject)B_mpi,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
374   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
375   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
376   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
377   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_STRUMPACK() is not implemented yet");
378   PetscFunctionReturn(0);
379 }
380 
381 #undef __FUNCT__
382 #define __FUNCT__ "MatFactorInfo_STRUMPACK"
383 static PetscErrorCode MatFactorInfo_STRUMPACK(Mat A,PetscViewer viewer)
384 {
385   PetscErrorCode  ierr;
386 
387   PetscFunctionBegin;
388   /* check if matrix is strumpack type */
389   if (A->ops->solve != MatSolve_STRUMPACK) PetscFunctionReturn(0);
390   ierr = PetscViewerASCIIPrintf(viewer,"STRUMPACK sparse solver!\n");CHKERRQ(ierr);
391   PetscFunctionReturn(0);
392 }
393 
394 #undef __FUNCT__
395 #define __FUNCT__ "MatView_STRUMPACK"
396 static PetscErrorCode MatView_STRUMPACK(Mat A,PetscViewer viewer)
397 {
398   PetscErrorCode    ierr;
399   PetscBool         iascii;
400   PetscViewerFormat format;
401 
402   PetscFunctionBegin;
403   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
404   if (iascii) {
405     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
406     if (format == PETSC_VIEWER_ASCII_INFO) {
407       ierr = MatFactorInfo_STRUMPACK(A,viewer);CHKERRQ(ierr);
408     }
409   }
410   PetscFunctionReturn(0);
411 }
412 
413 #undef __FUNCT__
414 #define __FUNCT__ "MatLUFactorNumeric_STRUMPACK"
415 static PetscErrorCode MatLUFactorNumeric_STRUMPACK(Mat F,Mat A,const MatFactorInfo *info)
416 {
417   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
418   STRUMPACK_RETURN_CODE  sp_err;
419   Mat_SeqAIJ             *A_d,*A_o;
420   Mat_MPIAIJ             *mat;
421   PetscErrorCode         ierr;
422   PetscInt               M=A->rmap->N,m=A->rmap->n;
423   PetscBool              flg;
424 
425   PetscFunctionBegin;
426   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
427   if (flg) { /* A is MATMPIAIJ */
428     mat = (Mat_MPIAIJ*)A->data;
429     A_d = (Mat_SeqAIJ*)(mat->A)->data;
430     A_o = (Mat_SeqAIJ*)(mat->B)->data;
431     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));
432   } else { /* A is MATSEQAIJ */
433     A_d = (Mat_SeqAIJ*)A->data;
434     PetscStackCall("STRUMPACK_set_csr_matrix",STRUMPACK_set_csr_matrix(*S,&M,A_d->i,A_d->j,A_d->a,0));
435   }
436 
437   /* Reorder and Factor the matrix. */
438   /* TODO figure out how to avoid reorder if the matrix values changed, but the pattern remains the same. */
439   PetscStackCall("STRUMPACK_reorder",sp_err = STRUMPACK_reorder(*S));
440   PetscStackCall("STRUMPACK_factor",sp_err = STRUMPACK_factor(*S));
441   switch (sp_err) {
442   case STRUMPACK_SUCCESS: break;
443   case STRUMPACK_MATRIX_NOT_SET:   { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); break; }
444   case STRUMPACK_REORDERING_ERROR: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); break; }
445   default:                           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: factorization failed");
446   }
447   PetscFunctionReturn(0);
448 }
449 
450 #undef __FUNCT__
451 #define __FUNCT__ "MatLUFactorSymbolic_STRUMPACK"
452 static PetscErrorCode MatLUFactorSymbolic_STRUMPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
453 {
454   PetscFunctionBegin;
455   F->ops->lufactornumeric = MatLUFactorNumeric_STRUMPACK;
456   F->ops->solve           = MatSolve_STRUMPACK;
457   F->ops->matsolve        = MatMatSolve_STRUMPACK;
458   PetscFunctionReturn(0);
459 }
460 
461 #undef __FUNCT__
462 #define __FUNCT__ "MatFactorGetSolverPackage_aij_strumpack"
463 static PetscErrorCode MatFactorGetSolverPackage_aij_strumpack(Mat A,const MatSolverPackage *type)
464 {
465   PetscFunctionBegin;
466   *type = MATSOLVERSTRUMPACK;
467   PetscFunctionReturn(0);
468 }
469 
470 /*MC
471   MATSOLVERSSTRUMPACK = "strumpack" - A solver package providing a direct sparse solver (PCLU)
472   and a preconditioner (PCILU) using low-rank compression via the external package STRUMPACK.
473 
474   Consult the STRUMPACK-sparse manual for more info.
475 
476   Use
477      ./configure --download-strumpack
478   to have PETSc installed with STRUMPACK
479 
480   Use
481     -pc_type lu -pc_factor_mat_solver_package strumpack
482   to use this as an exact (direct) solver, use
483     -pc_type ilu -pc_factor_mat_solver_package strumpack
484   to enable low-rank compression (i.e, use as a preconditioner).
485 
486   Works with AIJ matrices
487 
488   Options Database Keys:
489 + -mat_strumpack_verbose
490 . -mat_strumpack_hss_rel_tol <1e-2>         - Relative compression tolerance (None)
491 . -mat_strumpack_hss_abs_tol <1e-8>         - Absolute compression tolerance (None)
492 . -mat_strumpack_colperm <TRUE>             - Permute matrix to make diagonal nonzeros (None)
493 . -mat_strumpack_hss_min_front_size <1000>  - Minimum size of dense block for HSS compression (None)
494 . -mat_strumpack_hss_min_sep_size <256>     - Minimum size of separator for HSS compression (None)
495 . -mat_strumpack_max_rank                   - Maximum rank in HSS compression, when using pctype ilu (None)
496 . -mat_strumpack_reordering <METIS>         - Sparsity reducing matrix reordering (choose one of) NATURAL METIS PARMETIS SCOTCH PTSCOTCH RCM (None)
497 . -mat_strumpack_iterative_solver <DIRECT>  - Select iterative solver from STRUMPACK (choose one of) AUTO DIRECT REFINE PREC_GMRES GMRES PREC_BICGSTAB BICGSTAB (None)
498 
499  Level: beginner
500 
501 .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, PCFactorSetMatSolverPackage(), MatSolverPackage
502 M*/
503 #undef __FUNCT__
504 #define __FUNCT__ "MatGetFactor_aij_strumpack"
505 static PetscErrorCode MatGetFactor_aij_strumpack(Mat A,MatFactorType ftype,Mat *F)
506 {
507   Mat                           B;
508   PetscErrorCode                ierr;
509   PetscInt                      M=A->rmap->N,N=A->cmap->N;
510   PetscBool                     verb,flg,set;
511   PetscReal                     ctol;
512   PetscInt                      hssminsize,max_rank;
513   STRUMPACK_SparseSolver        *S;
514   STRUMPACK_INTERFACE           iface;
515   STRUMPACK_REORDERING_STRATEGY ndcurrent,ndvalue;
516   STRUMPACK_KRYLOV_SOLVER       itcurrent,itsolver;
517     const STRUMPACK_PRECISION table[2][2][2] =
518     {{{STRUMPACK_FLOATCOMPLEX_64, STRUMPACK_DOUBLECOMPLEX_64},
519       {STRUMPACK_FLOAT_64,        STRUMPACK_DOUBLE_64}},
520      {{STRUMPACK_FLOATCOMPLEX,    STRUMPACK_DOUBLECOMPLEX},
521       {STRUMPACK_FLOAT,           STRUMPACK_DOUBLE}}};
522   const STRUMPACK_PRECISION prec =
523     table[(sizeof(PetscInt)==8)?0:1]
524     [(PETSC_SCALAR==PETSC_COMPLEX)?0:1]
525     [(PETSC_REAL==PETSC_FLOAT)?0:1];
526 
527   PetscFunctionBegin;
528   /* Create the factorization matrix */
529   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
530   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,M,N);CHKERRQ(ierr);
531   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
532   ierr = MatSeqAIJSetPreallocation(B,0,NULL);
533   ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr);
534   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) {
535     B->ops->lufactorsymbolic  = MatLUFactorSymbolic_STRUMPACK;
536     B->ops->ilufactorsymbolic = MatLUFactorSymbolic_STRUMPACK;
537   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
538   B->ops->view        = MatView_STRUMPACK;
539   B->ops->destroy     = MatDestroy_STRUMPACK;
540   B->ops->getdiagonal = MatGetDiagonal_STRUMPACK;
541   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_aij_strumpack);CHKERRQ(ierr);
542   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetReordering_C",MatSTRUMPACKSetReordering_STRUMPACK);CHKERRQ(ierr);
543   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetColPerm_C",MatSTRUMPACKSetColPerm_STRUMPACK);CHKERRQ(ierr);
544   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSRelTol_C",MatSTRUMPACKSetHSSRelTol_STRUMPACK);CHKERRQ(ierr);
545   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSAbsTol_C",MatSTRUMPACKSetHSSAbsTol_STRUMPACK);CHKERRQ(ierr);
546   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSMinFrontSize_C",MatSTRUMPACKSetHSSMinFrontSize_STRUMPACK);CHKERRQ(ierr);
547   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSMinSepSize_C",MatSTRUMPACKSetHSSMinSepSize_STRUMPACK);CHKERRQ(ierr);
548   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSMaxRank_C",MatSTRUMPACKSetHSSMaxRank_STRUMPACK);CHKERRQ(ierr);
549   B->factortype = ftype;
550   ierr     = PetscNewLog(B,&S);CHKERRQ(ierr);
551   B->spptr = S;
552 
553   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);
554   iface = flg ? STRUMPACK_MT : STRUMPACK_MPI_DIST;
555 
556   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"STRUMPACK Options","Mat");CHKERRQ(ierr);
557 
558   verb = PetscLogPrintInfo ? PETSC_TRUE : PETSC_FALSE;
559   ierr = PetscOptionsBool("-mat_strumpack_verbose","Print STRUMPACK information","None",verb,&verb,NULL);CHKERRQ(ierr);
560 
561   PetscStackCall("STRUMPACK_init",STRUMPACK_init(S,PetscObjectComm((PetscObject)A),prec,iface,0,NULL,verb));
562 
563   PetscStackCall("STRUMPACK_hss_rel_tol",ctol = (PetscReal)STRUMPACK_hss_rel_tol(*S));
564   ierr = PetscOptionsReal("-mat_strumpack_hss_rel_tol","Relative compression tolerance","None",ctol,&ctol,&set);CHKERRQ(ierr);
565   if (set) PetscStackCall("STRUMPACK_set_hss_rel_tol",STRUMPACK_set_hss_rel_tol(*S,(double)ctol));
566 
567   PetscStackCall("STRUMPACK_hss_abs_tol",ctol = (PetscReal)STRUMPACK_hss_abs_tol(*S));
568   ierr = PetscOptionsReal("-mat_strumpack_hss_abs_tol","Absolute compression tolerance","None",ctol,&ctol,&set);CHKERRQ(ierr);
569   if (set) PetscStackCall("STRUMPACK_set_hss_abs_tol",STRUMPACK_set_hss_abs_tol(*S,(double)ctol));
570 
571   PetscStackCall("STRUMPACK_mc64job",flg = (STRUMPACK_mc64job(*S) == 0) ? PETSC_FALSE : PETSC_TRUE);
572   ierr = PetscOptionsBool("-mat_strumpack_colperm","Find a col perm to get nonzero diagonal","None",flg,&flg,&set);CHKERRQ(ierr);
573   if (set) PetscStackCall("STRUMPACK_set_mc64job",STRUMPACK_set_mc64job(*S,flg ? 5 : 0));
574 
575   PetscStackCall("STRUMPACK_HSS_min_front_size",hssminsize = (PetscInt)STRUMPACK_HSS_min_front_size(*S));
576   ierr = PetscOptionsInt("-mat_strumpack_hss_min_front_size","Minimum size of dense block for HSS compression","None",hssminsize,&hssminsize,&set);CHKERRQ(ierr);
577   if (set) PetscStackCall("STRUMPACK_set_HSS_min_front_size",STRUMPACK_set_HSS_min_front_size(*S,(int)hssminsize));
578 
579   PetscStackCall("STRUMPACK_HSS_min_sep_size",hssminsize = (PetscInt)STRUMPACK_HSS_min_sep_size(*S));
580   ierr = PetscOptionsInt("-mat_strumpack_hss_min_sep_size","Minimum size of separator for HSS compression","None",hssminsize,&hssminsize,&set);CHKERRQ(ierr);
581   if (set) PetscStackCall("STRUMPACK_set_HSS_min_sep_size",STRUMPACK_set_HSS_min_sep_size(*S,(int)hssminsize));
582 
583   PetscStackCall("STRUMPACK_max_rank",max_rank = (PetscInt)STRUMPACK_max_rank(*S));
584   ierr = PetscOptionsInt("-mat_strumpack_max_rank","Maximum rank in HSS approximation","None",max_rank,&max_rank,&set);CHKERRQ(ierr);
585   if (set) PetscStackCall("STRUMPACK_set_max_rank",STRUMPACK_set_max_rank(*S,(int)max_rank));
586 
587   const char *const STRUMPACKNDTypes[] = {"NATURAL","METIS","PARMETIS","SCOTCH","PTSCOTCH","RCM","STRUMPACKNDTypes","",0};
588   PetscStackCall("STRUMPACK_reordering_method",ndcurrent = STRUMPACK_reordering_method(*S));
589   PetscOptionsEnum("-mat_strumpack_reordering","Sparsity reducing matrix reordering","None",STRUMPACKNDTypes,(PetscEnum)ndcurrent,(PetscEnum*)&ndvalue,&set);CHKERRQ(ierr);
590   if (set) PetscStackCall("STRUMPACK_set_reordering_method",STRUMPACK_set_reordering_method(*S,ndvalue));
591 
592   if (ftype == MAT_FACTOR_ILU) {
593     /* When enabling HSS compression, the STRUMPACK solver becomes an incomplete              */
594     /* (or approximate) LU factorization.                                                     */
595     PetscStackCall("STRUMPACK_enable_HSS",STRUMPACK_enable_HSS(*S));
596   }
597 
598   /* Disable the outer iterative solver from STRUMPACK.                                       */
599   /* When STRUMPACK is used as a direct solver, it will by default do iterative refinement.   */
600   /* When STRUMPACK is used as an approximate factorization preconditioner (by enabling       */
601   /* low-rank compression), it will use it's own preconditioned GMRES. Here we can disable    */
602   /* the outer iterative solver, as PETSc uses STRUMPACK from within a KSP.                   */
603   PetscStackCall("STRUMPACK_set_Krylov_solver", STRUMPACK_set_Krylov_solver(*S, STRUMPACK_DIRECT));
604 
605   const char *const SolverTypes[] = {"AUTO","NONE","REFINE","PREC_GMRES","GMRES","PREC_BICGSTAB","BICGSTAB","SolverTypes","",0};
606   PetscStackCall("STRUMPACK_Krylov_solver",itcurrent = STRUMPACK_Krylov_solver(*S));
607   PetscOptionsEnum("-mat_strumpack_iterative_solver","Select iterative solver from STRUMPACK","None",SolverTypes,(PetscEnum)itcurrent,(PetscEnum*)&itsolver,&set);CHKERRQ(ierr);
608   if (set) PetscStackCall("STRUMPACK_set_Krylov_solver",STRUMPACK_set_Krylov_solver(*S,itsolver));
609 
610   PetscOptionsEnd();
611 
612   *F = B;
613   PetscFunctionReturn(0);
614 }
615 
616 #undef __FUNCT__
617 #define __FUNCT__ "MatSolverPackageRegister_STRUMPACK"
618 PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_STRUMPACK(void)
619 {
620   PetscErrorCode ierr;
621 
622   PetscFunctionBegin;
623   ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr);
624   ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr);
625   ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_ILU,MatGetFactor_aij_strumpack);CHKERRQ(ierr);
626   ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_ILU,MatGetFactor_aij_strumpack);CHKERRQ(ierr);
627   PetscFunctionReturn(0);
628 }
629