xref: /petsc/src/mat/impls/aij/seq/klu/klu.c (revision 69e15a41f829c371e4f5c4ef2075fdc7e71d4503)
1 
2 /*
3    Provides an interface to the KLUv1.2 sparse solver
4 
5    When build with PETSC_USE_64BIT_INDICES this will use SuiteSparse_long as the
6    integer type in KLU, otherwise it will use int. This means
7    all integers in this file are simply declared as PetscInt. Also it means
8    that KLU SuiteSparse_long version MUST be built with 64 bit integers when used.
9 
10 */
11 #include <../src/mat/impls/aij/seq/aij.h>
12 
13 #if defined(PETSC_USE_64BIT_INDICES)
14 #define klu_K_defaults                klu_l_defaults
15 #define klu_K_analyze                 klu_l_analyze
16 #define klu_K_analyze_given           klu_l_analyze_given
17 #define klu_K_free_symbolic           klu_l_free_symbolic
18 #define klu_K_free_numeric            klu_l_free_numeric
19 #define klu_K_common                  klu_l_common
20 #define klu_K_symbolic                klu_l_symbolic
21 #define klu_K_numeric                 klu_l_numeric
22 #if defined(PETSC_USE_COMPLEX)
23 #define klu_K_factor                  klu_zl_factor
24 #define klu_K_solve                   klu_zl_solve
25 #define klu_K_tsolve                  klu_zl_tsolve
26 #define klu_K_refactor                klu_zl_refactor
27 #define klu_K_sort                    klu_zl_sort
28 #define klu_K_flops                   klu_zl_flops
29 #define klu_K_rgrowth                 klu_zl_rgrowth
30 #define klu_K_condest                 klu_zl_condest
31 #define klu_K_rcond                   klu_zl_rcond
32 #define klu_K_scale                   klu_zl_scale
33 #else
34 #define klu_K_factor                  klu_l_factor
35 #define klu_K_solve                   klu_l_solve
36 #define klu_K_tsolve                  klu_l_tsolve
37 #define klu_K_refactor                klu_l_refactor
38 #define klu_K_sort                    klu_l_sort
39 #define klu_K_flops                   klu_l_flops
40 #define klu_K_rgrowth                 klu_l_rgrowth
41 #define klu_K_condest                 klu_l_condest
42 #define klu_K_rcond                   klu_l_rcond
43 #define klu_K_scale                   klu_l_scale
44 #endif
45 #else
46 #define klu_K_defaults                klu_defaults
47 #define klu_K_analyze                 klu_analyze
48 #define klu_K_analyze_given           klu_analyze_given
49 #define klu_K_free_symbolic           klu_free_symbolic
50 #define klu_K_free_numeric            klu_free_numeric
51 #define klu_K_common                  klu_common
52 #define klu_K_symbolic                klu_symbolic
53 #define klu_K_numeric                 klu_numeric
54 #if defined(PETSC_USE_COMPLEX)
55 #define klu_K_factor                  klu_z_factor
56 #define klu_K_solve                   klu_z_solve
57 #define klu_K_tsolve                  klu_z_tsolve
58 #define klu_K_refactor                klu_z_refactor
59 #define klu_K_sort                    klu_z_sort
60 #define klu_K_flops                   klu_z_flops
61 #define klu_K_rgrowth                 klu_z_rgrowth
62 #define klu_K_condest                 klu_z_condest
63 #define klu_K_rcond                   klu_z_rcond
64 #define klu_K_scale                   klu_z_scale
65 #else
66 #define klu_K_factor                  klu_factor
67 #define klu_K_solve                   klu_solve
68 #define klu_K_tsolve                  klu_tsolve
69 #define klu_K_refactor                klu_refactor
70 #define klu_K_sort                    klu_sort
71 #define klu_K_flops                   klu_flops
72 #define klu_K_rgrowth                 klu_rgrowth
73 #define klu_K_condest                 klu_condest
74 #define klu_K_rcond                   klu_rcond
75 #define klu_K_scale                   klu_scale
76 #endif
77 #endif
78 
79 
80 #define SuiteSparse_long long long
81 #define SuiteSparse_long_max LONG_LONG_MAX
82 #define SuiteSparse_long_id "%lld"
83 
84 EXTERN_C_BEGIN
85 #include <klu.h>
86 EXTERN_C_END
87 
88 static const char *KluOrderingTypes[] = {"AMD","COLAMD","PETSC"};
89 static const char *scale[] ={"NONE","SUM","MAX"};
90 
91 typedef struct {
92   klu_K_common   Common;
93   klu_K_symbolic *Symbolic;
94   klu_K_numeric  *Numeric;
95   PetscInt     *perm_c,*perm_r;
96   MatStructure flg;
97   PetscBool    PetscMatOrdering;
98 
99   /* Flag to clean up KLU objects during Destroy */
100   PetscBool CleanUpKLU;
101 } Mat_KLU;
102 
103 #undef __FUNCT__
104 #define __FUNCT__ "MatDestroy_KLU"
105 static PetscErrorCode MatDestroy_KLU(Mat A)
106 {
107   PetscErrorCode ierr;
108   Mat_KLU    *lu=(Mat_KLU*)A->spptr;
109 
110   PetscFunctionBegin;
111   if (lu && lu->CleanUpKLU) {
112     klu_K_free_symbolic(&lu->Symbolic,&lu->Common);
113     klu_K_free_numeric(&lu->Numeric,&lu->Common);
114     ierr = PetscFree(lu->perm_r);CHKERRQ(ierr);
115     ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
116   }
117   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
118   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
119   PetscFunctionReturn(0);
120 }
121 
122 #undef __FUNCT__
123 #define __FUNCT__ "MatSolveTranspose_KLU"
124 static PetscErrorCode MatSolveTranspose_KLU(Mat A,Vec b,Vec x)
125 {
126   Mat_KLU       *lu = (Mat_KLU*)A->spptr;
127   PetscScalar    *xa;
128   PetscErrorCode ierr;
129   PetscInt       status;
130 
131   PetscFunctionBegin;
132   /* KLU uses a column major format, solve Ax = b by klu_*_solve */
133   /* ----------------------------------*/
134   ierr = VecCopy(b,x); /* klu_solve stores the solution in rhs */
135   ierr = VecGetArray(x,&xa);
136   status = klu_K_solve(lu->Symbolic,lu->Numeric,A->rmap->n,1,(PetscReal*)xa,&lu->Common);
137   if (status != 1) {
138     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Solve failed");
139   }
140   ierr = VecRestoreArray(x,&xa);CHKERRQ(ierr);
141   PetscFunctionReturn(0);
142 }
143 
144 #undef __FUNCT__
145 #define __FUNCT__ "MatSolve_KLU"
146 static PetscErrorCode MatSolve_KLU(Mat A,Vec b,Vec x)
147 {
148   Mat_KLU       *lu = (Mat_KLU*)A->spptr;
149   PetscScalar    *xa;
150   PetscErrorCode ierr;
151   PetscInt       status;
152 
153   PetscFunctionBegin;
154   /* KLU uses a column major format, solve A^Tx = b by klu_*_tsolve */
155   /* ----------------------------------*/
156   ierr = VecCopy(b,x); /* klu_solve stores the solution in rhs */
157   ierr = VecGetArray(x,&xa);
158 #if defined(PETSC_USE_COMPLEX)
159   PetscInt conj_solve=1;
160   status = klu_K_tsolve(lu->Symbolic,lu->Numeric,A->rmap->n,1,(PetscReal*)xa,conj_solve,&lu->Common); /* conjugate solve */
161 #else
162   status = klu_K_tsolve(lu->Symbolic,lu->Numeric,A->rmap->n,1,xa,&lu->Common);
163 #endif
164   if (status != 1) {
165     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Solve failed");
166   }
167   ierr = VecRestoreArray(x,&xa);CHKERRQ(ierr);
168   PetscFunctionReturn(0);
169 }
170 
171 #undef __FUNCT__
172 #define __FUNCT__ "MatLUFactorNumeric_KLU"
173 static PetscErrorCode MatLUFactorNumeric_KLU(Mat F,Mat A,const MatFactorInfo *info)
174 {
175   Mat_KLU        *lu = (Mat_KLU*)(F)->spptr;
176   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
177   PetscInt       *ai = a->i,*aj=a->j;
178   PetscScalar    *av = a->a;
179 
180   PetscFunctionBegin;
181   /* numeric factorization of A' */
182   /* ----------------------------*/
183 
184   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric) {
185     klu_K_free_numeric(&lu->Numeric,&lu->Common);
186   }
187   lu->Numeric = klu_K_factor(ai,aj,(PetscReal*)av,lu->Symbolic,&lu->Common);
188   if(!lu->Numeric) {
189     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Numeric factorization failed");
190   }
191 
192   lu->flg                = SAME_NONZERO_PATTERN;
193   lu->CleanUpKLU         = PETSC_TRUE;
194   F->ops->solve          = MatSolve_KLU;
195   F->ops->solvetranspose = MatSolveTranspose_KLU;
196   PetscFunctionReturn(0);
197 }
198 
199 #undef __FUNCT__
200 #define __FUNCT__ "MatLUFactorSymbolic_KLU"
201 static PetscErrorCode MatLUFactorSymbolic_KLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
202 {
203   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
204   Mat_KLU       *lu = (Mat_KLU*)(F->spptr);
205   PetscErrorCode ierr;
206   PetscInt       i,*ai = a->i,*aj = a->j,m=A->rmap->n,n=A->cmap->n;
207   const PetscInt *ra,*ca;
208 
209   PetscFunctionBegin;
210   if (lu->PetscMatOrdering) {
211     ierr = ISGetIndices(r,&ra);CHKERRQ(ierr);
212     ierr = ISGetIndices(c,&ca);CHKERRQ(ierr);
213     ierr = PetscMalloc2(m,&lu->perm_c,n,&lu->perm_r);CHKERRQ(ierr);
214     /* we cannot simply memcpy on 64 bit archs */
215     for (i = 0; i < m; i++) lu->perm_r[i] = ra[i];
216     for (i = 0; i < n; i++) lu->perm_c[i] = ca[i];
217     ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr);
218     ierr = ISRestoreIndices(c,&ca);CHKERRQ(ierr);
219   }
220 
221   /* symbolic factorization of A' */
222   /* ---------------------------------------------------------------------- */
223   if (lu->PetscMatOrdering) { /* use Petsc ordering */
224     lu->Symbolic = klu_K_analyze_given(n,ai,aj,lu->perm_c,lu->perm_r,&lu->Common);
225   } else { /* use klu internal ordering */
226     lu->Symbolic = klu_K_analyze(n,ai,aj,&lu->Common);
227   }
228   if (!lu->Symbolic) {
229     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Symbolic Factorization failed");
230   }
231 
232   lu->flg                   = DIFFERENT_NONZERO_PATTERN;
233   lu->CleanUpKLU            = PETSC_TRUE;
234   (F)->ops->lufactornumeric = MatLUFactorNumeric_KLU;
235   PetscFunctionReturn(0);
236 }
237 
238 #undef __FUNCT__
239 #define __FUNCT__ "MatFactorInfo_KLU"
240 static PetscErrorCode MatFactorInfo_KLU(Mat A,PetscViewer viewer)
241 {
242   Mat_KLU       *lu= (Mat_KLU*)A->spptr;
243   klu_K_numeric *Numeric=(klu_K_numeric*)lu->Numeric;
244   PetscErrorCode ierr;
245 
246   PetscFunctionBegin;
247   /* check if matrix is KLU type */
248   if (A->ops->solve != MatSolve_KLU) PetscFunctionReturn(0);
249 
250   ierr = PetscViewerASCIIPrintf(viewer,"KLU stats:\n");CHKERRQ(ierr);
251   ierr = PetscViewerASCIIPrintf(viewer,"  Number of diagonal blocks: %d\n",Numeric->nblocks);
252   ierr = PetscViewerASCIIPrintf(viewer,"  Total nonzeros=%d\n",Numeric->lnz+Numeric->unz);CHKERRQ(ierr);
253 
254   ierr = PetscViewerASCIIPrintf(viewer,"KLU runtime parameters:\n");CHKERRQ(ierr);
255 
256   /* Control parameters used by numeric factorization */
257   ierr = PetscViewerASCIIPrintf(viewer,"  Partial pivoting tolerance: %g\n",lu->Common.tol);CHKERRQ(ierr);
258   /* BTF preordering */
259   ierr = PetscViewerASCIIPrintf(viewer,"  BTF preordering enabled: %d\n",lu->Common.btf);CHKERRQ(ierr);
260   /* mat ordering */
261   if (!lu->PetscMatOrdering) {
262     ierr = PetscViewerASCIIPrintf(viewer,"  Ordering: %s (not using the PETSc ordering)\n",KluOrderingTypes[(int)lu->Common.ordering]);CHKERRQ(ierr);
263   } else {
264     ierr = PetscViewerASCIIPrintf(viewer,"  Using PETSc ordering\n");CHKERRQ(ierr);
265   }
266   /* matrix row scaling */
267   ierr = PetscViewerASCIIPrintf(viewer, "  Matrix row scaling: %s\n",scale[(int)lu->Common.scale]);CHKERRQ(ierr);
268   PetscFunctionReturn(0);
269 }
270 
271 #undef __FUNCT__
272 #define __FUNCT__ "MatView_KLU"
273 static PetscErrorCode MatView_KLU(Mat A,PetscViewer viewer)
274 {
275   PetscErrorCode    ierr;
276   PetscBool         iascii;
277   PetscViewerFormat format;
278 
279   PetscFunctionBegin;
280   ierr = MatView_SeqAIJ(A,viewer);CHKERRQ(ierr);
281 
282   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
283   if (iascii) {
284     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
285     if (format == PETSC_VIEWER_ASCII_INFO) {
286       ierr = MatFactorInfo_KLU(A,viewer);CHKERRQ(ierr);
287     }
288   }
289   PetscFunctionReturn(0);
290 }
291 
292 #undef __FUNCT__
293 #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_klu"
294 PetscErrorCode MatFactorGetSolverPackage_seqaij_klu(Mat A,const MatSolverPackage *type)
295 {
296   PetscFunctionBegin;
297   *type = MATSOLVERKLU;
298   PetscFunctionReturn(0);
299 }
300 
301 
302 /*MC
303   MATSOLVERKLU = "klu" - A matrix type providing direct solvers (LU) for sequential matrices
304   via the external package KLU.
305 
306   ./configure --download-SuiteSparse --SuiteSparse-enable-klu to install PETSc to use KLU
307 
308   Consult KLU documentation for more information on the options database keys below.
309 
310   Options Database Keys:
311 + -mat_klu_pivot_tol <0.001>                  - Partial pivoting tolerance
312 . -mat_klu_use_btf <1>                        - Use BTF preordering
313 . -mat_klu_ordering <AMD>                     - KLU reordering scheme to reduce fill-in (choose one of) AMD COLAMD PETSC
314 - -mat_klu_row_scale <NONE>                   - Matrix row scaling (choose one of) NONE SUM MAX
315 
316    Level: beginner
317 
318 .seealso: PCLU, MATSOLVERUMFPACK, MATSOLVERCHOLMOD, PCFactorSetMatSolverPackage(), MatSolverPackage
319 M*/
320 
321 #undef __FUNCT__
322 #define __FUNCT__ "MatGetFactor_seqaij_klu"
323 PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_klu(Mat A,MatFactorType ftype,Mat *F)
324 {
325   Mat            B;
326   Mat_KLU       *lu;
327   PetscErrorCode ierr;
328   PetscInt       m=A->rmap->n,n=A->cmap->n,idx,status;
329   PetscBool      flg;
330 
331   PetscFunctionBegin;
332   /* Create the factorization matrix F */
333   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
334   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
335   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
336   ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr);
337   ierr = PetscNewLog(B,&lu);CHKERRQ(ierr);
338 
339   B->spptr                 = lu;
340   B->ops->lufactorsymbolic = MatLUFactorSymbolic_KLU;
341   B->ops->destroy          = MatDestroy_KLU;
342   B->ops->view             = MatView_KLU;
343 
344   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_seqaij_klu);CHKERRQ(ierr);
345 
346   B->factortype   = MAT_FACTOR_LU;
347   B->assembled    = PETSC_TRUE;           /* required by -ksp_view */
348   B->preallocated = PETSC_TRUE;
349 
350   /* initializations */
351   /* ------------------------------------------------*/
352   /* get the default control parameters */
353   status = klu_K_defaults(&lu->Common);
354   if(status <= 0) {
355     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Initialization failed");
356   }
357   lu->Common.scale = 0; /* No row scaling */
358 
359   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"KLU Options","Mat");CHKERRQ(ierr);
360   /* Partial pivoting tolerance */
361   ierr = PetscOptionsReal("-mat_klu_pivot_tol","Partial pivoting tolerance","None",lu->Common.tol,&lu->Common.tol,NULL);CHKERRQ(ierr);
362   /* BTF pre-ordering */
363   ierr = PetscOptionsInt("-mat_klu_use_btf","Enable BTF preordering","None",lu->Common.btf,&lu->Common.btf,NULL);CHKERRQ(ierr);
364   /* Matrix reordering */
365   ierr = PetscOptionsEList("-mat_klu_ordering","Internal ordering method","None",KluOrderingTypes,sizeof(KluOrderingTypes)/sizeof(KluOrderingTypes[0]),KluOrderingTypes[0],&idx,&flg);CHKERRQ(ierr);
366   if (flg) {
367     if ((int)idx == 2) lu->PetscMatOrdering = PETSC_TRUE;   /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Klu perm_c) */
368     else lu->Common.ordering = (int)idx;
369   }
370   /* Matrix row scaling */
371   ierr = PetscOptionsEList("-mat_klu_row_scale","Matrix row scaling","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr);
372   PetscOptionsEnd();
373   *F = B;
374   PetscFunctionReturn(0);
375 }
376