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