xref: /petsc/src/mat/impls/aij/mpi/strumpack/strumpack.c (revision 609bdbee21ea3be08735c64dbe00a9ab27759925)
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,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr);
31   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSRelCompTol_C",NULL);CHKERRQ(ierr);
32   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSMinSize_C",NULL);CHKERRQ(ierr);
33   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetColPerm_C",NULL);CHKERRQ(ierr);
34 
35   PetscFunctionReturn(0);
36 }
37 
38 static PetscErrorCode MatSTRUMPACKSetHSSRelCompTol_STRUMPACK(Mat F,PetscReal rctol)
39 {
40   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
41 
42   PetscFunctionBegin;
43   PetscStackCall("STRUMPACK_set_rctol", STRUMPACK_set_rctol(*S,rctol));
44   PetscFunctionReturn(0);
45 }
46 
47 /*@
48   MatSTRUMPACKSetHSSRelCompTol - Set STRUMPACK relative tolerance for HSS compression
49    Logically Collective on Mat
50 
51    Input Parameters:
52 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
53 -  rctol - relative compression tolerance
54 
55   Options Database:
56 .   -mat_strumpack_rctol <rctol>
57 
58    Level: beginner
59 
60    References:
61 .      STRUMPACK manual
62 
63 .seealso: MatGetFactor()
64 @*/
65 PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat F,PetscReal rctol)
66 {
67   PetscErrorCode ierr;
68 
69   PetscFunctionBegin;
70   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
71   PetscValidLogicalCollectiveReal(F,rctol,2);
72   ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSRelCompTol_C",(Mat,PetscReal),(F,rctol));CHKERRQ(ierr);
73   PetscFunctionReturn(0);
74 }
75 
76 static PetscErrorCode MatSTRUMPACKSetHSSMinSize_STRUMPACK(Mat F,PetscInt hssminsize)
77 {
78   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
79 
80   PetscFunctionBegin;
81   PetscStackCall("STRUMPACK_set_minimum_HSS_size", STRUMPACK_set_minimum_HSS_size(*S,hssminsize));
82   PetscFunctionReturn(0);
83 }
84 
85 /*@
86   MatSTRUMPACKSetHSSMinSize - Set STRUMPACK minimum dense matrix size for low-rank approximation
87    Logically Collective on Mat
88 
89    Input Parameters:
90 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
91 -  hssminsize - minimum dense matrix size for low-rank approximation
92 
93   Options Database:
94 .   -mat_strumpack_hssminsize <hssminsize>
95 
96    Level: beginner
97 
98    References:
99 .      STRUMPACK manual
100 
101 .seealso: MatGetFactor()
102 @*/
103 PetscErrorCode MatSTRUMPACKSetHSSMinSize(Mat F,PetscInt hssminsize)
104 {
105   PetscErrorCode ierr;
106 
107   PetscFunctionBegin;
108   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
109   PetscValidLogicalCollectiveInt(F,hssminsize,2);
110   ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSMinSize_C",(Mat,PetscInt),(F,hssminsize));CHKERRQ(ierr);
111   PetscFunctionReturn(0);
112 }
113 
114 static PetscErrorCode MatSTRUMPACKSetColPerm_STRUMPACK(Mat F,PetscBool cperm)
115 {
116   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
117 
118   PetscFunctionBegin;
119   PetscStackCall("STRUMPACK_set_mc64job",STRUMPACK_set_mc64job(*S,cperm ? 5 : 0));
120   PetscFunctionReturn(0);
121 }
122 
123 /*@
124   MatSTRUMPACKSetColPerm - Set whether STRUMPACK should try to permute the columns of the matrix in order to get a nonzero diagonal
125    Logically Collective on Mat
126 
127    Input Parameters:
128 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface
129 -  cperm - whether or not to permute (internally) the columns of the matrix
130 
131   Options Database:
132 .   -mat_strumpack_colperm <cperm>
133 
134    Level: beginner
135 
136    References:
137 .      STRUMPACK manual
138 
139 .seealso: MatGetFactor()
140 @*/
141 PetscErrorCode MatSTRUMPACKSetColPerm(Mat F,PetscBool cperm)
142 {
143   PetscErrorCode ierr;
144 
145   PetscFunctionBegin;
146   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
147   PetscValidLogicalCollectiveBool(F,cperm,2);
148   ierr = PetscTryMethod(F,"MatSTRUMPACKSetColPerm_C",(Mat,PetscBool),(F,cperm));CHKERRQ(ierr);
149   PetscFunctionReturn(0);
150 }
151 
152 static PetscErrorCode MatSolve_STRUMPACK(Mat A,Vec b_mpi,Vec x)
153 {
154   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)A->spptr;
155   STRUMPACK_RETURN_CODE  sp_err;
156   PetscErrorCode         ierr;
157   const PetscScalar      *bptr;
158   PetscScalar            *xptr;
159 
160   PetscFunctionBegin;
161   ierr = VecGetArray(x,&xptr);CHKERRQ(ierr);
162   ierr = VecGetArrayRead(b_mpi,&bptr);CHKERRQ(ierr);
163 
164   PetscStackCall("STRUMPACK_solve",sp_err = STRUMPACK_solve(*S,(PetscScalar*)bptr,xptr,0));
165   switch (sp_err) {
166   case STRUMPACK_SUCCESS: break;
167   case STRUMPACK_MATRIX_NOT_SET:   { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); break; }
168   case STRUMPACK_REORDERING_ERROR: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); break; }
169   default:                           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: solve failed");
170   }
171   ierr = VecRestoreArray(x,&xptr);CHKERRQ(ierr);
172   ierr = VecRestoreArrayRead(b_mpi,&bptr);CHKERRQ(ierr);
173   PetscFunctionReturn(0);
174 }
175 
176 static PetscErrorCode MatMatSolve_STRUMPACK(Mat A,Mat B_mpi,Mat X)
177 {
178   PetscErrorCode   ierr;
179   PetscBool        flg;
180 
181   PetscFunctionBegin;
182   ierr = PetscObjectTypeCompareAny((PetscObject)B_mpi,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
183   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
184   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
185   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
186   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_STRUMPACK() is not implemented yet");
187   PetscFunctionReturn(0);
188 }
189 
190 static PetscErrorCode MatFactorInfo_STRUMPACK(Mat A,PetscViewer viewer)
191 {
192   PetscErrorCode  ierr;
193 
194   PetscFunctionBegin;
195   /* check if matrix is strumpack type */
196   if (A->ops->solve != MatSolve_STRUMPACK) PetscFunctionReturn(0);
197   ierr = PetscViewerASCIIPrintf(viewer,"STRUMPACK sparse solver!\n");CHKERRQ(ierr);
198   PetscFunctionReturn(0);
199 }
200 
201 static PetscErrorCode MatView_STRUMPACK(Mat A,PetscViewer viewer)
202 {
203   PetscErrorCode    ierr;
204   PetscBool         iascii;
205   PetscViewerFormat format;
206 
207   PetscFunctionBegin;
208   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
209   if (iascii) {
210     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
211     if (format == PETSC_VIEWER_ASCII_INFO) {
212       ierr = MatFactorInfo_STRUMPACK(A,viewer);CHKERRQ(ierr);
213     }
214   }
215   PetscFunctionReturn(0);
216 }
217 
218 static PetscErrorCode MatLUFactorNumeric_STRUMPACK(Mat F,Mat A,const MatFactorInfo *info)
219 {
220   STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr;
221   STRUMPACK_RETURN_CODE  sp_err;
222   Mat_SeqAIJ             *A_d,*A_o;
223   Mat_MPIAIJ             *mat;
224   PetscErrorCode         ierr;
225   PetscInt               M=A->rmap->N,m=A->rmap->n;
226   PetscBool              flg;
227 
228   PetscFunctionBegin;
229   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
230   if (flg) { /* A is MATMPIAIJ */
231     mat = (Mat_MPIAIJ*)A->data;
232     A_d = (Mat_SeqAIJ*)(mat->A)->data;
233     A_o = (Mat_SeqAIJ*)(mat->B)->data;
234     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));
235   } else { /* A is MATSEQAIJ */
236     A_d = (Mat_SeqAIJ*)A->data;
237     PetscStackCall("STRUMPACK_set_csr_matrix",STRUMPACK_set_csr_matrix(*S,&M,A_d->i,A_d->j,A_d->a,0));
238   }
239 
240   /* Reorder and Factor the matrix. */
241   /* TODO figure out how to avoid reorder if the matrix values changed, but the pattern remains the same. */
242   PetscStackCall("STRUMPACK_reorder",sp_err = STRUMPACK_reorder(*S));
243   PetscStackCall("STRUMPACK_factor",sp_err = STRUMPACK_factor(*S));
244   switch (sp_err) {
245   case STRUMPACK_SUCCESS: break;
246   case STRUMPACK_MATRIX_NOT_SET:   { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); break; }
247   case STRUMPACK_REORDERING_ERROR: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); break; }
248   default:                           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: factorization failed");
249   }
250   PetscFunctionReturn(0);
251 }
252 
253 static PetscErrorCode MatLUFactorSymbolic_STRUMPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
254 {
255   PetscFunctionBegin;
256   F->ops->lufactornumeric = MatLUFactorNumeric_STRUMPACK;
257   F->ops->solve           = MatSolve_STRUMPACK;
258   F->ops->matsolve        = MatMatSolve_STRUMPACK;
259   PetscFunctionReturn(0);
260 }
261 
262 static PetscErrorCode MatFactorGetSolverPackage_aij_strumpack(Mat A,const MatSolverPackage *type)
263 {
264   PetscFunctionBegin;
265   *type = MATSOLVERSTRUMPACK;
266   PetscFunctionReturn(0);
267 }
268 
269 /*MC
270   MATSOLVERSSTRUMPACK = "strumpack" - A solver package providing a direct sparse solver (PCLU)
271   and a preconditioner (PCILU) using low-rank compression via the external package STRUMPACK.
272 
273   Consult the STRUMPACK-sparse manual for more info.
274 
275   Use
276      ./configure --download-strumpack
277   to have PETSc installed with STRUMPACK
278 
279   Use
280     -pc_type lu -pc_factor_mat_solver_package strumpack
281   to use this as an exact (direct) solver, use
282     -pc_type ilu -pc_factor_mat_solver_package strumpack
283   to enable low-rank compression (i.e, use as a preconditioner).
284 
285   Works with AIJ matrices
286 
287   Options Database Keys:
288 + -mat_strumpack_rctol <1e-4>           - Relative compression tolerance (None)
289 . -mat_strumpack_hssminsize <512>       - Minimum size of dense block for HSS compression (None)
290 . -mat_strumpack_colperm <TRUE>         - Permute matrix to make diagonal nonzeros (None)
291 
292  Level: beginner
293 
294 .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, PCFactorSetMatSolverPackage(), MatSolverPackage
295 M*/
296 static PetscErrorCode MatGetFactor_aij_strumpack(Mat A,MatFactorType ftype,Mat *F)
297 {
298   Mat                    B;
299   STRUMPACK_SparseSolver *S;
300   PetscErrorCode         ierr;
301   PetscInt               M=A->rmap->N,N=A->cmap->N;
302   STRUMPACK_INTERFACE    iface;
303   PetscBool              verb,flg,set;
304   PetscReal              rctol;
305   PetscInt               hssminsize;
306   int                    argc;
307   char                   **args,*copts,*pname;
308   size_t                 len;
309   const STRUMPACK_PRECISION table[2][2][2] = {{{STRUMPACK_FLOATCOMPLEX_64,STRUMPACK_DOUBLECOMPLEX_64},
310                                                {STRUMPACK_FLOAT_64,STRUMPACK_DOUBLE_64}},
311                                               {{STRUMPACK_FLOATCOMPLEX,STRUMPACK_DOUBLECOMPLEX},
312                                                {STRUMPACK_FLOAT,STRUMPACK_DOUBLE}}};
313   const STRUMPACK_PRECISION prec = table[(sizeof(PetscInt)==8)?0:1][(PETSC_SCALAR==PETSC_COMPLEX)?0:1][(PETSC_REAL==PETSC_FLOAT)?0:1];
314 
315   PetscFunctionBegin;
316   /* Create the factorization matrix */
317   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
318   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,M,N);CHKERRQ(ierr);
319   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
320   ierr = MatSeqAIJSetPreallocation(B,0,NULL);
321   ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr);
322   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) {
323     B->ops->lufactorsymbolic  = MatLUFactorSymbolic_STRUMPACK;
324     B->ops->ilufactorsymbolic = MatLUFactorSymbolic_STRUMPACK;
325   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
326   B->ops->view        = MatView_STRUMPACK;
327   B->ops->destroy     = MatDestroy_STRUMPACK;
328   B->ops->getdiagonal = MatGetDiagonal_STRUMPACK;
329   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_aij_strumpack);CHKERRQ(ierr);
330   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSRelCompTol_C",MatSTRUMPACKSetHSSRelCompTol_STRUMPACK);CHKERRQ(ierr);
331   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSMinSize_C",MatSTRUMPACKSetHSSMinSize_STRUMPACK);CHKERRQ(ierr);
332   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetColPerm_C",MatSTRUMPACKSetColPerm_STRUMPACK);CHKERRQ(ierr);
333   B->factortype = ftype;
334   ierr     = PetscNewLog(B,&S);CHKERRQ(ierr);
335   B->spptr = S;
336 
337   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);
338   iface = flg ? STRUMPACK_MT : STRUMPACK_MPI_DIST;
339 
340   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"STRUMPACK Options","Mat");CHKERRQ(ierr);
341 
342   verb = PetscLogPrintInfo ? PETSC_TRUE : PETSC_FALSE;
343   ierr = PetscOptionsBool("-mat_strumpack_verbose","Print STRUMPACK information","None",verb,&verb,NULL);CHKERRQ(ierr);
344 
345   ierr = PetscOptionsGetAll(NULL,&copts);CHKERRQ(ierr);
346   ierr = PetscStrlen(copts,&len);CHKERRQ(ierr);
347   len += PETSC_MAX_PATH_LEN+1;
348   ierr = PetscMalloc1(len,&pname);CHKERRQ(ierr);
349   /* first string is assumed to be the program name, so add program name to options string */
350   ierr = PetscGetProgramName(pname,len);CHKERRQ(ierr);
351   ierr = PetscStrcat(pname," ");CHKERRQ(ierr);
352   ierr = PetscStrcat(pname,copts);CHKERRQ(ierr);
353   ierr = PetscStrToArray(pname,' ',&argc,&args);CHKERRQ(ierr);
354   ierr = PetscFree(copts);CHKERRQ(ierr);
355   ierr = PetscFree(pname);CHKERRQ(ierr);
356 
357   PetscStackCall("STRUMPACK_init",STRUMPACK_init(S,PetscObjectComm((PetscObject)A),prec,iface,argc,args,verb));
358 
359   PetscStackCall("STRUMPACK_get_rctol",rctol = (PetscReal)STRUMPACK_get_rctol(*S));
360   ierr = PetscOptionsReal("-mat_strumpack_rctol","Relative compression tolerance","None",rctol,&rctol,&set);CHKERRQ(ierr);
361   if (set) PetscStackCall("STRUMPACK_set_rctol",STRUMPACK_set_rctol(*S,(double)rctol));
362 
363   PetscStackCall("STRUMPACK_get_mc64job",flg = (STRUMPACK_get_mc64job(*S) == 0) ? PETSC_FALSE : PETSC_TRUE);
364   ierr = PetscOptionsBool("-mat_strumpack_colperm","Find a col perm to get nonzero diagonal","None",flg,&flg,&set);CHKERRQ(ierr);
365   if (set) PetscStackCall("STRUMPACK_set_mc64job",STRUMPACK_set_mc64job(*S,flg ? 5 : 0));
366 
367   PetscStackCall("STRUMPACK_get_minimum_HSS_size",hssminsize = (PetscInt)STRUMPACK_get_minimum_HSS_size(*S));
368   ierr = PetscOptionsInt("-mat_strumpack_hssminsize","Minimum size of dense block for HSS compression","None",hssminsize,&hssminsize,&set);CHKERRQ(ierr);
369   if (set) PetscStackCall("STRUMPACK_set_minimum_HSS_size",STRUMPACK_set_minimum_HSS_size(*S,(int)hssminsize));
370 
371   PetscOptionsEnd();
372 
373   if (ftype == MAT_FACTOR_ILU) {
374     /* When enabling HSS compression, the STRUMPACK solver becomes an incomplete                */
375     /* (or approximate) LU factorization.                                                       */
376     PetscStackCall("STRUMPACK_use_HSS",STRUMPACK_use_HSS(*S,1));
377     /* Disable the outer iterative solver from STRUMPACK.                                       */
378     /* When STRUMPACK is used as a direct solver, it will by default do iterative refinement.   */
379     /* When STRUMPACK is used with as an approximate factorization preconditioner (by enabling  */
380     /* low-rank compression), it will use it's own GMRES. Here we can disable the               */
381     /* outer iterative solver, as PETSc uses STRUMPACK from within a KSP.                       */
382     PetscStackCall("STRUMPACK_set_Krylov_solver", STRUMPACK_set_Krylov_solver(*S, STRUMPACK_DIRECT));
383   }
384 
385   /* Allow the user to set or overwrite the above options from the command line                 */
386   PetscStackCall("STRUMPACK_set_from_options",STRUMPACK_set_from_options(*S));
387   ierr = PetscStrToArrayDestroy(argc,args);CHKERRQ(ierr);
388 
389   *F = B;
390   PetscFunctionReturn(0);
391 }
392 
393 PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_STRUMPACK(void)
394 {
395   PetscErrorCode ierr;
396 
397   PetscFunctionBegin;
398   ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr);
399   ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr);
400   ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_ILU,MatGetFactor_aij_strumpack);CHKERRQ(ierr);
401   ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_ILU,MatGetFactor_aij_strumpack);CHKERRQ(ierr);
402   PetscFunctionReturn(0);
403 }
404