xref: /petsc/src/mat/impls/aij/seq/klu/klu.c (revision cf1aed2ce99d23e50336629af3ca8cf096900abb)
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 = PetscFree2(lu->perm_r,lu->perm_c);CHKERRQ(ierr);
115   }
116   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
117   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
118   PetscFunctionReturn(0);
119 }
120 
121 #undef __FUNCT__
122 #define __FUNCT__ "MatSolveTranspose_KLU"
123 static PetscErrorCode MatSolveTranspose_KLU(Mat A,Vec b,Vec x)
124 {
125   Mat_KLU       *lu = (Mat_KLU*)A->spptr;
126   PetscScalar    *xa;
127   PetscErrorCode ierr;
128   PetscInt       status;
129 
130   PetscFunctionBegin;
131   /* KLU uses a column major format, solve Ax = b by klu_*_solve */
132   /* ----------------------------------*/
133   ierr = VecCopy(b,x); /* klu_solve stores the solution in rhs */
134   ierr = VecGetArray(x,&xa);
135   status = klu_K_solve(lu->Symbolic,lu->Numeric,A->rmap->n,1,(PetscReal*)xa,&lu->Common);
136   if (status != 1) {
137     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Solve failed");
138   }
139   ierr = VecRestoreArray(x,&xa);CHKERRQ(ierr);
140   PetscFunctionReturn(0);
141 }
142 
143 #undef __FUNCT__
144 #define __FUNCT__ "MatSolve_KLU"
145 static PetscErrorCode MatSolve_KLU(Mat A,Vec b,Vec x)
146 {
147   Mat_KLU       *lu = (Mat_KLU*)A->spptr;
148   PetscScalar    *xa;
149   PetscErrorCode ierr;
150   PetscInt       status;
151 
152   PetscFunctionBegin;
153   /* KLU uses a column major format, solve A^Tx = b by klu_*_tsolve */
154   /* ----------------------------------*/
155   ierr = VecCopy(b,x); /* klu_solve stores the solution in rhs */
156   ierr = VecGetArray(x,&xa);
157 #if defined(PETSC_USE_COMPLEX)
158   PetscInt conj_solve=1;
159   status = klu_K_tsolve(lu->Symbolic,lu->Numeric,A->rmap->n,1,(PetscReal*)xa,conj_solve,&lu->Common); /* conjugate solve */
160 #else
161   status = klu_K_tsolve(lu->Symbolic,lu->Numeric,A->rmap->n,1,xa,&lu->Common);
162 #endif
163   if (status != 1) {
164     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Solve failed");
165   }
166   ierr = VecRestoreArray(x,&xa);CHKERRQ(ierr);
167   PetscFunctionReturn(0);
168 }
169 
170 #undef __FUNCT__
171 #define __FUNCT__ "MatLUFactorNumeric_KLU"
172 static PetscErrorCode MatLUFactorNumeric_KLU(Mat F,Mat A,const MatFactorInfo *info)
173 {
174   Mat_KLU        *lu = (Mat_KLU*)(F)->spptr;
175   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
176   PetscInt       *ai = a->i,*aj=a->j;
177   PetscScalar    *av = a->a;
178 
179   PetscFunctionBegin;
180   /* numeric factorization of A' */
181   /* ----------------------------*/
182 
183   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric) {
184     klu_K_free_numeric(&lu->Numeric,&lu->Common);
185   }
186   lu->Numeric = klu_K_factor(ai,aj,(PetscReal*)av,lu->Symbolic,&lu->Common);
187   if(!lu->Numeric) {
188     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Numeric factorization failed");
189   }
190 
191   lu->flg                = SAME_NONZERO_PATTERN;
192   lu->CleanUpKLU         = PETSC_TRUE;
193   F->ops->solve          = MatSolve_KLU;
194   F->ops->solvetranspose = MatSolveTranspose_KLU;
195   PetscFunctionReturn(0);
196 }
197 
198 #undef __FUNCT__
199 #define __FUNCT__ "MatLUFactorSymbolic_KLU"
200 static PetscErrorCode MatLUFactorSymbolic_KLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
201 {
202   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
203   Mat_KLU       *lu = (Mat_KLU*)(F->spptr);
204   PetscErrorCode ierr;
205   PetscInt       i,*ai = a->i,*aj = a->j,m=A->rmap->n,n=A->cmap->n;
206   const PetscInt *ra,*ca;
207 
208   PetscFunctionBegin;
209   if (lu->PetscMatOrdering) {
210     ierr = ISGetIndices(r,&ra);CHKERRQ(ierr);
211     ierr = ISGetIndices(c,&ca);CHKERRQ(ierr);
212     ierr = PetscMalloc2(m,&lu->perm_r,n,&lu->perm_c);CHKERRQ(ierr);
213     /* we cannot simply memcpy on 64 bit archs */
214     for (i = 0; i < m; i++) lu->perm_r[i] = ra[i];
215     for (i = 0; i < n; i++) lu->perm_c[i] = ca[i];
216     ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr);
217     ierr = ISRestoreIndices(c,&ca);CHKERRQ(ierr);
218   }
219 
220   /* symbolic factorization of A' */
221   /* ---------------------------------------------------------------------- */
222   if (lu->PetscMatOrdering) { /* use Petsc ordering */
223     lu->Symbolic = klu_K_analyze_given(n,ai,aj,lu->perm_c,lu->perm_r,&lu->Common);
224   } else { /* use klu internal ordering */
225     lu->Symbolic = klu_K_analyze(n,ai,aj,&lu->Common);
226   }
227   if (!lu->Symbolic) {
228     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Symbolic Factorization failed");
229   }
230 
231   lu->flg                   = DIFFERENT_NONZERO_PATTERN;
232   lu->CleanUpKLU            = PETSC_TRUE;
233   (F)->ops->lufactornumeric = MatLUFactorNumeric_KLU;
234   PetscFunctionReturn(0);
235 }
236 
237 #undef __FUNCT__
238 #define __FUNCT__ "MatFactorInfo_KLU"
239 static PetscErrorCode MatFactorInfo_KLU(Mat A,PetscViewer viewer)
240 {
241   Mat_KLU       *lu= (Mat_KLU*)A->spptr;
242   klu_K_numeric *Numeric=(klu_K_numeric*)lu->Numeric;
243   PetscErrorCode ierr;
244 
245   PetscFunctionBegin;
246   /* check if matrix is KLU type */
247   if (A->ops->solve != MatSolve_KLU) PetscFunctionReturn(0);
248 
249   ierr = PetscViewerASCIIPrintf(viewer,"KLU stats:\n");CHKERRQ(ierr);
250   ierr = PetscViewerASCIIPrintf(viewer,"  Number of diagonal blocks: %d\n",Numeric->nblocks);
251   ierr = PetscViewerASCIIPrintf(viewer,"  Total nonzeros=%d\n",Numeric->lnz+Numeric->unz);CHKERRQ(ierr);
252 
253   ierr = PetscViewerASCIIPrintf(viewer,"KLU runtime parameters:\n");CHKERRQ(ierr);
254 
255   /* Control parameters used by numeric factorization */
256   ierr = PetscViewerASCIIPrintf(viewer,"  Partial pivoting tolerance: %g\n",lu->Common.tol);CHKERRQ(ierr);
257   /* BTF preordering */
258   ierr = PetscViewerASCIIPrintf(viewer,"  BTF preordering enabled: %d\n",lu->Common.btf);CHKERRQ(ierr);
259   /* mat ordering */
260   if (!lu->PetscMatOrdering) {
261     ierr = PetscViewerASCIIPrintf(viewer,"  Ordering: %s (not using the PETSc ordering)\n",KluOrderingTypes[(int)lu->Common.ordering]);CHKERRQ(ierr);
262   } else {
263     ierr = PetscViewerASCIIPrintf(viewer,"  Using PETSc ordering\n");CHKERRQ(ierr);
264   }
265   /* matrix row scaling */
266   ierr = PetscViewerASCIIPrintf(viewer, "  Matrix row scaling: %s\n",scale[(int)lu->Common.scale]);CHKERRQ(ierr);
267   PetscFunctionReturn(0);
268 }
269 
270 #undef __FUNCT__
271 #define __FUNCT__ "MatView_KLU"
272 static PetscErrorCode MatView_KLU(Mat A,PetscViewer viewer)
273 {
274   PetscErrorCode    ierr;
275   PetscBool         iascii;
276   PetscViewerFormat format;
277 
278   PetscFunctionBegin;
279   ierr = MatView_SeqAIJ(A,viewer);CHKERRQ(ierr);
280 
281   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
282   if (iascii) {
283     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
284     if (format == PETSC_VIEWER_ASCII_INFO) {
285       ierr = MatFactorInfo_KLU(A,viewer);CHKERRQ(ierr);
286     }
287   }
288   PetscFunctionReturn(0);
289 }
290 
291 #undef __FUNCT__
292 #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_klu"
293 PetscErrorCode MatFactorGetSolverPackage_seqaij_klu(Mat A,const MatSolverPackage *type)
294 {
295   PetscFunctionBegin;
296   *type = MATSOLVERKLU;
297   PetscFunctionReturn(0);
298 }
299 
300 
301 /*MC
302   MATSOLVERKLU = "klu" - A matrix type providing direct solvers (LU) for sequential matrices
303   via the external package KLU.
304 
305   ./configure --download-suitesparse to install PETSc to use KLU
306 
307   Consult KLU documentation for more information on the options database keys below.
308 
309   Options Database Keys:
310 + -mat_klu_pivot_tol <0.001>                  - Partial pivoting tolerance
311 . -mat_klu_use_btf <1>                        - Use BTF preordering
312 . -mat_klu_ordering <AMD>                     - KLU reordering scheme to reduce fill-in (choose one of) AMD COLAMD PETSC
313 - -mat_klu_row_scale <NONE>                   - Matrix row scaling (choose one of) NONE SUM MAX
314 
315    Level: beginner
316 
317 .seealso: PCLU, MATSOLVERUMFPACK, MATSOLVERCHOLMOD, PCFactorSetMatSolverPackage(), MatSolverPackage
318 M*/
319 
320 #undef __FUNCT__
321 #define __FUNCT__ "MatGetFactor_seqaij_klu"
322 PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_klu(Mat A,MatFactorType ftype,Mat *F)
323 {
324   Mat            B;
325   Mat_KLU       *lu;
326   PetscErrorCode ierr;
327   PetscInt       m=A->rmap->n,n=A->cmap->n,idx,status;
328   PetscBool      flg;
329 
330   PetscFunctionBegin;
331   /* Create the factorization matrix F */
332   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
333   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
334   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
335   ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr);
336   ierr = PetscNewLog(B,&lu);CHKERRQ(ierr);
337 
338   B->spptr                 = lu;
339   B->ops->lufactorsymbolic = MatLUFactorSymbolic_KLU;
340   B->ops->destroy          = MatDestroy_KLU;
341   B->ops->view             = MatView_KLU;
342 
343   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_seqaij_klu);CHKERRQ(ierr);
344 
345   B->factortype   = MAT_FACTOR_LU;
346   B->assembled    = PETSC_TRUE;           /* required by -ksp_view */
347   B->preallocated = PETSC_TRUE;
348 
349   /* initializations */
350   /* ------------------------------------------------*/
351   /* get the default control parameters */
352   status = klu_K_defaults(&lu->Common);
353   if(status <= 0) {
354     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Initialization failed");
355   }
356   lu->Common.scale = 0; /* No row scaling */
357 
358   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"KLU Options","Mat");CHKERRQ(ierr);
359   /* Partial pivoting tolerance */
360   ierr = PetscOptionsReal("-mat_klu_pivot_tol","Partial pivoting tolerance","None",lu->Common.tol,&lu->Common.tol,NULL);CHKERRQ(ierr);
361   /* BTF pre-ordering */
362   ierr = PetscOptionsInt("-mat_klu_use_btf","Enable BTF preordering","None",lu->Common.btf,&lu->Common.btf,NULL);CHKERRQ(ierr);
363   /* Matrix reordering */
364   ierr = PetscOptionsEList("-mat_klu_ordering","Internal ordering method","None",KluOrderingTypes,sizeof(KluOrderingTypes)/sizeof(KluOrderingTypes[0]),KluOrderingTypes[0],&idx,&flg);CHKERRQ(ierr);
365   if (flg) {
366     if ((int)idx == 2) lu->PetscMatOrdering = PETSC_TRUE;   /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Klu perm_c) */
367     else lu->Common.ordering = (int)idx;
368   }
369   /* Matrix row scaling */
370   ierr = PetscOptionsEList("-mat_klu_row_scale","Matrix row scaling","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr);
371   PetscOptionsEnd();
372   *F = B;
373   PetscFunctionReturn(0);
374 }
375