xref: /petsc/src/mat/impls/aij/seq/klu/klu.c (revision 7d5fd1e4d9337468ad3f05b65b7facdcd2dfd2a4)
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(a,b,c,d)           klu_l_analyze((SuiteSparse_long)a,(SuiteSparse_long*)b,(SuiteSparse_long*)c,d)
16 #define klu_K_analyze_given(a,b,c,d,e,f) klu_l_analyze_given((SuiteSparse_long)a,(SuiteSparse_long*)b,(SuiteSparse_long*)c,(SuiteSparse_long*)d,(SuiteSparse_long*)e,f)
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(a,b,c,d,e)       klu_zl_factor((SuiteSparse_long*)a,(SuiteSparse_long*)b,c,d,e);
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(a,b,c,d,e)       klu_l_factor((SuiteSparse_long*)a,(SuiteSparse_long*)b,c,d,e);
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 EXTERN_C_BEGIN
80 #include <klu.h>
81 EXTERN_C_END
82 
83 static const char *KluOrderingTypes[] = {"AMD","COLAMD"};
84 static const char *scale[] ={"NONE","SUM","MAX"};
85 
86 typedef struct {
87   klu_K_common   Common;
88   klu_K_symbolic *Symbolic;
89   klu_K_numeric  *Numeric;
90   PetscInt       *perm_c,*perm_r;
91   MatStructure   flg;
92   PetscBool      PetscMatOrdering;
93   PetscBool      CleanUpKLU;
94 } Mat_KLU;
95 
96 static PetscErrorCode MatDestroy_KLU(Mat A)
97 {
98   PetscErrorCode ierr;
99   Mat_KLU    *lu=(Mat_KLU*)A->data;
100 
101   PetscFunctionBegin;
102   if (lu->CleanUpKLU) {
103     klu_K_free_symbolic(&lu->Symbolic,&lu->Common);
104     klu_K_free_numeric(&lu->Numeric,&lu->Common);
105     ierr = PetscFree2(lu->perm_r,lu->perm_c);CHKERRQ(ierr);
106   }
107   ierr = PetscFree(A->data);CHKERRQ(ierr);
108   PetscFunctionReturn(0);
109 }
110 
111 static PetscErrorCode MatSolveTranspose_KLU(Mat A,Vec b,Vec x)
112 {
113   Mat_KLU       *lu = (Mat_KLU*)A->data;
114   PetscScalar    *xa;
115   PetscErrorCode ierr;
116   PetscInt       status;
117 
118   PetscFunctionBegin;
119   /* KLU uses a column major format, solve Ax = b by klu_*_solve */
120   /* ----------------------------------*/
121   ierr = VecCopy(b,x); /* klu_solve stores the solution in rhs */
122   ierr = VecGetArray(x,&xa);
123   status = klu_K_solve(lu->Symbolic,lu->Numeric,A->rmap->n,1,(PetscReal*)xa,&lu->Common);
124   if (status != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Solve failed");
125   ierr = VecRestoreArray(x,&xa);CHKERRQ(ierr);
126   PetscFunctionReturn(0);
127 }
128 
129 static PetscErrorCode MatSolve_KLU(Mat A,Vec b,Vec x)
130 {
131   Mat_KLU       *lu = (Mat_KLU*)A->data;
132   PetscScalar    *xa;
133   PetscErrorCode ierr;
134   PetscInt       status;
135 
136   PetscFunctionBegin;
137   /* KLU uses a column major format, solve A^Tx = b by klu_*_tsolve */
138   /* ----------------------------------*/
139   ierr = VecCopy(b,x); /* klu_solve stores the solution in rhs */
140   ierr = VecGetArray(x,&xa);
141 #if defined(PETSC_USE_COMPLEX)
142   PetscInt conj_solve=1;
143   status = klu_K_tsolve(lu->Symbolic,lu->Numeric,A->rmap->n,1,(PetscReal*)xa,conj_solve,&lu->Common); /* conjugate solve */
144 #else
145   status = klu_K_tsolve(lu->Symbolic,lu->Numeric,A->rmap->n,1,xa,&lu->Common);
146 #endif
147   if (status != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Solve failed");
148   ierr = VecRestoreArray(x,&xa);CHKERRQ(ierr);
149   PetscFunctionReturn(0);
150 }
151 
152 static PetscErrorCode MatLUFactorNumeric_KLU(Mat F,Mat A,const MatFactorInfo *info)
153 {
154   Mat_KLU        *lu = (Mat_KLU*)(F)->data;
155   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
156   PetscInt       *ai = a->i,*aj=a->j;
157   PetscScalar    *av = a->a;
158 
159   PetscFunctionBegin;
160   /* numeric factorization of A' */
161   /* ----------------------------*/
162 
163   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric) {
164     klu_K_free_numeric(&lu->Numeric,&lu->Common);
165   }
166   lu->Numeric = klu_K_factor(ai,aj,(PetscReal*)av,lu->Symbolic,&lu->Common);
167   if (!lu->Numeric) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Numeric factorization failed");
168 
169   lu->flg                = SAME_NONZERO_PATTERN;
170   lu->CleanUpKLU         = PETSC_TRUE;
171   F->ops->solve          = MatSolve_KLU;
172   F->ops->solvetranspose = MatSolveTranspose_KLU;
173   PetscFunctionReturn(0);
174 }
175 
176 static PetscErrorCode MatLUFactorSymbolic_KLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
177 {
178   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
179   Mat_KLU        *lu = (Mat_KLU*)(F->data);
180   PetscErrorCode ierr;
181   PetscInt       i,*ai = a->i,*aj = a->j,m=A->rmap->n,n=A->cmap->n;
182   const PetscInt *ra,*ca;
183 
184   PetscFunctionBegin;
185   if (lu->PetscMatOrdering) {
186     ierr = ISGetIndices(r,&ra);CHKERRQ(ierr);
187     ierr = ISGetIndices(c,&ca);CHKERRQ(ierr);
188     ierr = PetscMalloc2(m,&lu->perm_r,n,&lu->perm_c);CHKERRQ(ierr);
189     /* we cannot simply memcpy on 64 bit archs */
190     for (i = 0; i < m; i++) lu->perm_r[i] = ra[i];
191     for (i = 0; i < n; i++) lu->perm_c[i] = ca[i];
192     ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr);
193     ierr = ISRestoreIndices(c,&ca);CHKERRQ(ierr);
194   }
195 
196   /* symbolic factorization of A' */
197   /* ---------------------------------------------------------------------- */
198   if (r) {
199     lu->PetscMatOrdering = PETSC_TRUE;
200     lu->Symbolic = klu_K_analyze_given(n,ai,aj,lu->perm_c,lu->perm_r,&lu->Common);
201   } else { /* use klu internal ordering */
202     lu->Symbolic = klu_K_analyze(n,ai,aj,&lu->Common);
203   }
204   if (!lu->Symbolic) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Symbolic Factorization failed");
205 
206   lu->flg                   = DIFFERENT_NONZERO_PATTERN;
207   lu->CleanUpKLU            = PETSC_TRUE;
208   (F)->ops->lufactornumeric = MatLUFactorNumeric_KLU;
209   PetscFunctionReturn(0);
210 }
211 
212 static PetscErrorCode MatView_Info_KLU(Mat A,PetscViewer viewer)
213 {
214   Mat_KLU       *lu= (Mat_KLU*)A->data;
215   klu_K_numeric *Numeric=(klu_K_numeric*)lu->Numeric;
216   PetscErrorCode ierr;
217 
218   PetscFunctionBegin;
219   ierr = PetscViewerASCIIPrintf(viewer,"KLU stats:\n");CHKERRQ(ierr);
220   ierr = PetscViewerASCIIPrintf(viewer,"  Number of diagonal blocks: %d\n",Numeric->nblocks);
221   ierr = PetscViewerASCIIPrintf(viewer,"  Total nonzeros=%d\n",Numeric->lnz+Numeric->unz);CHKERRQ(ierr);
222   ierr = PetscViewerASCIIPrintf(viewer,"KLU runtime parameters:\n");CHKERRQ(ierr);
223   /* Control parameters used by numeric factorization */
224   ierr = PetscViewerASCIIPrintf(viewer,"  Partial pivoting tolerance: %g\n",lu->Common.tol);CHKERRQ(ierr);
225   /* BTF preordering */
226   ierr = PetscViewerASCIIPrintf(viewer,"  BTF preordering enabled: %d\n",lu->Common.btf);CHKERRQ(ierr);
227   /* mat ordering */
228   if (!lu->PetscMatOrdering) {
229     ierr = PetscViewerASCIIPrintf(viewer,"  Ordering: %s (not using the PETSc ordering)\n",KluOrderingTypes[(int)lu->Common.ordering]);CHKERRQ(ierr);
230   }
231   /* matrix row scaling */
232   ierr = PetscViewerASCIIPrintf(viewer, "  Matrix row scaling: %s\n",scale[(int)lu->Common.scale]);CHKERRQ(ierr);
233   PetscFunctionReturn(0);
234 }
235 
236 static PetscErrorCode MatView_KLU(Mat A,PetscViewer viewer)
237 {
238   PetscErrorCode    ierr;
239   PetscBool         iascii;
240   PetscViewerFormat format;
241 
242   PetscFunctionBegin;
243   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
244   if (iascii) {
245     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
246     if (format == PETSC_VIEWER_ASCII_INFO) {
247       ierr = MatView_Info_KLU(A,viewer);CHKERRQ(ierr);
248     }
249   }
250   PetscFunctionReturn(0);
251 }
252 
253 PetscErrorCode MatFactorGetSolverType_seqaij_klu(Mat A,MatSolverType *type)
254 {
255   PetscFunctionBegin;
256   *type = MATSOLVERKLU;
257   PetscFunctionReturn(0);
258 }
259 
260 /*MC
261   MATSOLVERKLU = "klu" - A matrix type providing direct solvers (LU) for sequential matrices
262   via the external package KLU.
263 
264   ./configure --download-suitesparse to install PETSc to use KLU
265 
266   Use -pc_type lu -pc_factor_mat_solver_type klu to use this direct solver
267 
268   Consult KLU documentation for more information on the options database keys below.
269 
270   Options Database Keys:
271 + -mat_klu_pivot_tol <0.001>                  - Partial pivoting tolerance
272 . -mat_klu_use_btf <1>                        - Use BTF preordering
273 . -mat_klu_ordering <AMD>                     - KLU reordering scheme to reduce fill-in (choose one of) AMD COLAMD PETSC
274 - -mat_klu_row_scale <NONE>                   - Matrix row scaling (choose one of) NONE SUM MAX
275 
276    Note: KLU is part of SuiteSparse http://faculty.cse.tamu.edu/davis/suitesparse.html
277 
278    Level: beginner
279 
280 .seealso: PCLU, MATSOLVERUMFPACK, MATSOLVERCHOLMOD, PCFactorSetMatSolverType(), MatSolverType
281 M*/
282 
283 PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_klu(Mat A,MatFactorType ftype,Mat *F)
284 {
285   Mat            B;
286   Mat_KLU       *lu;
287   PetscErrorCode ierr;
288   PetscInt       m=A->rmap->n,n=A->cmap->n,idx = 0,status;
289   PetscBool      flg;
290 
291   PetscFunctionBegin;
292   /* Create the factorization matrix F */
293   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
294   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
295   ierr = PetscStrallocpy("klu",&((PetscObject)B)->type_name);CHKERRQ(ierr);
296   ierr = MatSetUp(B);CHKERRQ(ierr);
297 
298   ierr = PetscNewLog(B,&lu);CHKERRQ(ierr);
299 
300   B->data                  = lu;
301   B->ops->getinfo          = MatGetInfo_External;
302   B->ops->lufactorsymbolic = MatLUFactorSymbolic_KLU;
303   B->ops->destroy          = MatDestroy_KLU;
304   B->ops->view             = MatView_KLU;
305 
306   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_seqaij_klu);CHKERRQ(ierr);
307 
308   B->factortype   = MAT_FACTOR_LU;
309   B->assembled    = PETSC_TRUE;           /* required by -ksp_view */
310   B->preallocated = PETSC_TRUE;
311 
312   ierr = PetscFree(B->solvertype);CHKERRQ(ierr);
313   ierr = PetscStrallocpy(MATSOLVERKLU,&B->solvertype);CHKERRQ(ierr);
314   B->canuseordering = PETSC_TRUE;
315   ierr = PetscStrallocpy(MATORDERINGEXTERNAL,(char**)&B->preferredordering[MAT_FACTOR_LU]);CHKERRQ(ierr);
316 
317   /* initializations */
318   /* ------------------------------------------------*/
319   /* get the default control parameters */
320   status = klu_K_defaults(&lu->Common);
321   if (status <= 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Initialization failed");
322 
323   lu->Common.scale = 0; /* No row scaling */
324 
325   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"KLU Options","Mat");CHKERRQ(ierr);
326   /* Partial pivoting tolerance */
327   ierr = PetscOptionsReal("-mat_klu_pivot_tol","Partial pivoting tolerance","None",lu->Common.tol,&lu->Common.tol,NULL);CHKERRQ(ierr);
328   /* BTF pre-ordering */
329   ierr = PetscOptionsInt("-mat_klu_use_btf","Enable BTF preordering","None",(PetscInt)lu->Common.btf,(PetscInt*)&lu->Common.btf,NULL);CHKERRQ(ierr);
330   /* Matrix reordering */
331   ierr = PetscOptionsEList("-mat_klu_ordering","Internal ordering method","None",KluOrderingTypes,sizeof(KluOrderingTypes)/sizeof(KluOrderingTypes[0]),KluOrderingTypes[0],&idx,&flg);CHKERRQ(ierr);
332   lu->Common.ordering = (int)idx;
333   /* Matrix row scaling */
334   ierr = PetscOptionsEList("-mat_klu_row_scale","Matrix row scaling","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr);
335   PetscOptionsEnd();
336   *F = B;
337   PetscFunctionReturn(0);
338 }
339